2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/math/functions.h"
44 #include "gromacs/mdlib/gmx_omp_nthreads.h"
45 #include "gromacs/mdlib/qmmm.h"
46 #include "gromacs/mdtypes/inputrec.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/topology/mtop_lookup.h"
49 #include "gromacs/topology/mtop_util.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/alignedallocator.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/smalloc.h"
55 #define ALMOST_ZERO 1e-30
57 t_mdatoms *init_mdatoms(FILE *fp, const gmx_mtop_t &mtop, const t_inputrec &ir)
62 md->nenergrp = mtop.groups.grps[egcENER].nr;
63 md->bVCMgrps = (mtop.groups.grps[egcVCM].nr > 1);
65 /* Determine the total system mass and perturbed atom counts */
66 double totalMassA = 0.0;
67 double totalMassB = 0.0;
69 gmx_mtop_atomloop_block_t aloop = gmx_mtop_atomloop_block_init(&mtop);
72 while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
74 totalMassA += nmol*atom->m;
75 totalMassB += nmol*atom->mB;
77 if (ir.efep != efepNO && PERTURBED(*atom))
80 if (atom->mB != atom->m)
82 md->nMassPerturbed += nmol;
84 if (atom->qB != atom->q)
86 md->nChargePerturbed += nmol;
88 if (atom->typeB != atom->type)
90 md->nTypePerturbed += nmol;
95 md->tmassA = totalMassA;
96 md->tmassB = totalMassB;
98 if (ir.efep != efepNO && fp)
101 "There are %d atoms and %d charges for free energy perturbation\n",
102 md->nPerturbed, md->nChargePerturbed);
105 md->havePartiallyFrozenAtoms = FALSE;
106 for (int g = 0; g < ir.opts.ngfrz; g++)
108 for (int d = YY; d < DIM; d++)
110 if (ir.opts.nFreeze[d] != ir.opts.nFreeze[XX])
112 md->havePartiallyFrozenAtoms = TRUE;
117 md->bOrires = gmx_mtop_ftype_count(&mtop, F_ORIRES);
122 void atoms2md(const gmx_mtop_t *mtop, const t_inputrec *ir,
123 int nindex, const int *index,
128 const t_grpopts *opts;
129 const gmx_groups_t *groups;
130 int nthreads gmx_unused;
132 bLJPME = EVDW_PME(ir->vdwtype);
136 groups = &mtop->groups;
138 /* nindex>=0 indicates DD where we use an index */
145 md->nr = mtop->natoms;
148 if (md->nr > md->nalloc)
150 md->nalloc = over_alloc_dd(md->nr);
152 if (md->nMassPerturbed)
154 srenew(md->massA, md->nalloc);
155 srenew(md->massB, md->nalloc);
157 srenew(md->massT, md->nalloc);
158 /* The SIMD version of the integrator needs this aligned and padded.
159 * The padding needs to be with zeros, which we set later below.
161 gmx::AlignedAllocationPolicy::free(md->invmass);
162 md->invmass = new(gmx::AlignedAllocationPolicy::malloc((md->nalloc + GMX_REAL_MAX_SIMD_WIDTH)*sizeof(*md->invmass)))real;
163 srenew(md->invMassPerDim, md->nalloc);
164 srenew(md->chargeA, md->nalloc);
165 srenew(md->typeA, md->nalloc);
168 srenew(md->chargeB, md->nalloc);
169 srenew(md->typeB, md->nalloc);
173 srenew(md->sqrt_c6A, md->nalloc);
174 srenew(md->sigmaA, md->nalloc);
175 srenew(md->sigma3A, md->nalloc);
178 srenew(md->sqrt_c6B, md->nalloc);
179 srenew(md->sigmaB, md->nalloc);
180 srenew(md->sigma3B, md->nalloc);
183 srenew(md->ptype, md->nalloc);
186 srenew(md->cTC, md->nalloc);
187 /* We always copy cTC with domain decomposition */
189 srenew(md->cENER, md->nalloc);
192 srenew(md->cACC, md->nalloc);
196 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
198 srenew(md->cFREEZE, md->nalloc);
202 srenew(md->cVCM, md->nalloc);
206 srenew(md->cORF, md->nalloc);
210 srenew(md->bPerturbed, md->nalloc);
213 /* Note that these user t_mdatoms array pointers are NULL
214 * when there is only one group present.
215 * Therefore, when adding code, the user should use something like:
216 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
218 if (mtop->groups.grpnr[egcUser1] != nullptr)
220 srenew(md->cU1, md->nalloc);
222 if (mtop->groups.grpnr[egcUser2] != nullptr)
224 srenew(md->cU2, md->nalloc);
229 srenew(md->bQM, md->nalloc);
235 // cppcheck-suppress unreadVariable
236 nthreads = gmx_omp_nthreads_get(emntDefault);
237 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
238 for (int i = 0; i < md->nr; i++)
246 if (index == nullptr)
254 const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
258 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
260 if (EI_ENERGY_MINIMIZATION(ir->eI))
262 /* Displacement is proportional to F, masses used for constraints */
266 else if (ir->eI == eiBD)
268 /* With BD the physical masses are irrelevant.
269 * To keep the code simple we use most of the normal MD code path
270 * for BD. Thus for constraining the masses should be proportional
271 * to the friction coefficient. We set the absolute value such that
272 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
273 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
274 * correct kinetic energy and temperature using the usual code path.
275 * Thus with BD v*dt will give the displacement and the reported
276 * temperature can signal bad integration (too large time step).
280 mA = 0.5*ir->bd_fric*ir->delta_t;
281 mB = 0.5*ir->bd_fric*ir->delta_t;
285 /* The friction coefficient is mass/tau_t */
286 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
288 mB = 0.5*atom.mB*fac;
296 if (md->nMassPerturbed)
306 md->invMassPerDim[i][XX] = 0;
307 md->invMassPerDim[i][YY] = 0;
308 md->invMassPerDim[i][ZZ] = 0;
310 else if (md->cFREEZE)
313 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
315 /* Set the mass of completely frozen particles to ALMOST_ZERO
316 * iso 0 to avoid div by zero in lincs or shake.
318 md->invmass[i] = ALMOST_ZERO;
322 /* Note: Partially frozen particles use the normal invmass.
323 * If such particles are constrained, the frozen dimensions
324 * should not be updated with the constrained coordinates.
326 md->invmass[i] = 1.0/mA;
328 for (int d = 0; d < DIM; d++)
330 md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0/mA);
335 md->invmass[i] = 1.0/mA;
336 for (int d = 0; d < DIM; d++)
338 md->invMassPerDim[i][d] = 1.0/mA;
342 md->chargeA[i] = atom.q;
343 md->typeA[i] = atom.type;
346 c6 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c6;
347 c12 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c12;
348 md->sqrt_c6A[i] = sqrt(c6);
349 if (c6 == 0.0 || c12 == 0)
355 md->sigmaA[i] = gmx::sixthroot(c12/c6);
357 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
361 md->bPerturbed[i] = PERTURBED(atom);
362 md->chargeB[i] = atom.qB;
363 md->typeB[i] = atom.typeB;
366 c6 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c6;
367 c12 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c12;
368 md->sqrt_c6B[i] = sqrt(c6);
369 if (c6 == 0.0 || c12 == 0)
375 md->sigmaB[i] = gmx::sixthroot(c12/c6);
377 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
380 md->ptype[i] = atom.ptype;
383 md->cTC[i] = groups->grpnr[egcTC][ag];
385 md->cENER[i] = ggrpnr(groups, egcENER, ag);
388 md->cACC[i] = groups->grpnr[egcACC][ag];
392 md->cVCM[i] = groups->grpnr[egcVCM][ag];
396 md->cORF[i] = ggrpnr(groups, egcORFIT, ag);
401 md->cU1[i] = groups->grpnr[egcUser1][ag];
405 md->cU2[i] = groups->grpnr[egcUser2][ag];
410 if (groups->grpnr[egcQMMM] == nullptr ||
411 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
421 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
426 /* Pad invmass with 0 so a SIMD MD update does not change v and x */
427 for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
434 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
435 * For free-energy runs, these should be updated using update_mdatoms().
437 md->tmass = md->tmassA;
441 void update_mdatoms(t_mdatoms *md, real lambda)
443 if (md->nMassPerturbed && lambda != md->lambda)
445 real L1 = 1 - lambda;
447 /* Update masses of perturbed atoms for the change in lambda */
448 // cppcheck-suppress unreadVariable
449 int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
450 #pragma omp parallel for num_threads(nthreads) schedule(static)
451 for (int i = 0; i < md->nr; i++)
453 if (md->bPerturbed[i])
455 md->massT[i] = L1*md->massA[i] + lambda*md->massB[i];
456 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
457 * and their invmass does not depend on lambda.
459 if (md->invmass[i] > 1.1*ALMOST_ZERO)
461 md->invmass[i] = 1.0/md->massT[i];
462 for (int d = 0; d < DIM; d++)
464 if (md->invMassPerDim[i][d] > 1.1*ALMOST_ZERO)
466 md->invMassPerDim[i][d] = md->invmass[i];
473 /* Update the system mass for the change in lambda */
474 md->tmass = L1*md->tmassA + lambda*md->tmassB;