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/exceptions.h"
52 #include "gromacs/utility/smalloc.h"
54 #define ALMOST_ZERO 1e-30
56 t_mdatoms *init_mdatoms(FILE *fp, const gmx_mtop_t *mtop, gmx_bool bFreeEnergy)
62 gmx_mtop_atomloop_all_t aloop;
66 md->nenergrp = mtop->groups.grps[egcENER].nr;
71 aloop = gmx_mtop_atomloop_all_init(mtop);
72 while (gmx_mtop_atomloop_all_next(aloop, &a, &atom))
74 if (ggrpnr(&mtop->groups, egcVCM, a) > 0)
79 if (bFreeEnergy && PERTURBED(*atom))
82 if (atom->mB != atom->m)
86 if (atom->qB != atom->q)
88 md->nChargePerturbed++;
90 if (atom->typeB != atom->type)
103 if (bFreeEnergy && fp)
106 "There are %d atoms and %d charges for free energy perturbation\n",
107 md->nPerturbed, md->nChargePerturbed);
110 md->bOrires = gmx_mtop_ftype_count(mtop, F_ORIRES);
115 void atoms2md(const gmx_mtop_t *mtop, const t_inputrec *ir,
116 int nindex, const int *index,
121 const t_grpopts *opts;
122 const gmx_groups_t *groups;
123 int nthreads gmx_unused;
125 bLJPME = EVDW_PME(ir->vdwtype);
129 groups = &mtop->groups;
131 /* nindex>=0 indicates DD where we use an index */
138 md->nr = mtop->natoms;
141 if (md->nr > md->nalloc)
143 md->nalloc = over_alloc_dd(md->nr);
145 if (md->nMassPerturbed)
147 srenew(md->massA, md->nalloc);
148 srenew(md->massB, md->nalloc);
150 srenew(md->massT, md->nalloc);
151 srenew(md->invmass, md->nalloc);
152 srenew(md->invMassPerDim, md->nalloc);
153 srenew(md->chargeA, md->nalloc);
154 srenew(md->typeA, md->nalloc);
157 srenew(md->chargeB, md->nalloc);
158 srenew(md->typeB, md->nalloc);
162 srenew(md->sqrt_c6A, md->nalloc);
163 srenew(md->sigmaA, md->nalloc);
164 srenew(md->sigma3A, md->nalloc);
167 srenew(md->sqrt_c6B, md->nalloc);
168 srenew(md->sigmaB, md->nalloc);
169 srenew(md->sigma3B, md->nalloc);
172 srenew(md->ptype, md->nalloc);
175 srenew(md->cTC, md->nalloc);
176 /* We always copy cTC with domain decomposition */
178 srenew(md->cENER, md->nalloc);
181 srenew(md->cACC, md->nalloc);
185 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
187 srenew(md->cFREEZE, md->nalloc);
191 srenew(md->cVCM, md->nalloc);
195 srenew(md->cORF, md->nalloc);
199 srenew(md->bPerturbed, md->nalloc);
202 /* Note that these user t_mdatoms array pointers are NULL
203 * when there is only one group present.
204 * Therefore, when adding code, the user should use something like:
205 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
207 if (mtop->groups.grpnr[egcUser1] != nullptr)
209 srenew(md->cU1, md->nalloc);
211 if (mtop->groups.grpnr[egcUser2] != nullptr)
213 srenew(md->cU2, md->nalloc);
218 srenew(md->bQM, md->nalloc);
224 // cppcheck-suppress unreadVariable
225 nthreads = gmx_omp_nthreads_get(emntDefault);
226 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
227 for (int i = 0; i < md->nr; i++)
235 if (index == nullptr)
243 const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
247 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
249 if (EI_ENERGY_MINIMIZATION(ir->eI))
251 /* Displacement is proportional to F, masses used for constraints */
255 else if (ir->eI == eiBD)
257 /* With BD the physical masses are irrelevant.
258 * To keep the code simple we use most of the normal MD code path
259 * for BD. Thus for constraining the masses should be proportional
260 * to the friction coefficient. We set the absolute value such that
261 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
262 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
263 * correct kinetic energy and temperature using the usual code path.
264 * Thus with BD v*dt will give the displacement and the reported
265 * temperature can signal bad integration (too large time step).
269 mA = 0.5*ir->bd_fric*ir->delta_t;
270 mB = 0.5*ir->bd_fric*ir->delta_t;
274 /* The friction coefficient is mass/tau_t */
275 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
277 mB = 0.5*atom.mB*fac;
285 if (md->nMassPerturbed)
295 md->invMassPerDim[i][XX] = 0;
296 md->invMassPerDim[i][YY] = 0;
297 md->invMassPerDim[i][ZZ] = 0;
299 else if (md->cFREEZE)
302 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
304 /* Set the mass of completely frozen particles to ALMOST_ZERO
305 * iso 0 to avoid div by zero in lincs or shake.
307 md->invmass[i] = ALMOST_ZERO;
311 /* Note: Partially frozen particles use the normal invmass.
312 * If such particles are constrained, the frozen dimensions
313 * should not be updated with the constrained coordinates.
315 md->invmass[i] = 1.0/mA;
317 for (int d = 0; d < DIM; d++)
319 md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0/mA);
324 md->invmass[i] = 1.0/mA;
325 for (int d = 0; d < DIM; d++)
327 md->invMassPerDim[i][d] = 1.0/mA;
331 md->chargeA[i] = atom.q;
332 md->typeA[i] = atom.type;
335 c6 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c6;
336 c12 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c12;
337 md->sqrt_c6A[i] = sqrt(c6);
338 if (c6 == 0.0 || c12 == 0)
344 md->sigmaA[i] = gmx::sixthroot(c12/c6);
346 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
350 md->bPerturbed[i] = PERTURBED(atom);
351 md->chargeB[i] = atom.qB;
352 md->typeB[i] = atom.typeB;
355 c6 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c6;
356 c12 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c12;
357 md->sqrt_c6B[i] = sqrt(c6);
358 if (c6 == 0.0 || c12 == 0)
364 md->sigmaB[i] = gmx::sixthroot(c12/c6);
366 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
369 md->ptype[i] = atom.ptype;
372 md->cTC[i] = groups->grpnr[egcTC][ag];
375 (groups->grpnr[egcENER] ? groups->grpnr[egcENER][ag] : 0);
378 md->cACC[i] = groups->grpnr[egcACC][ag];
382 md->cVCM[i] = groups->grpnr[egcVCM][ag];
386 md->cORF[i] = groups->grpnr[egcORFIT][ag];
391 md->cU1[i] = groups->grpnr[egcUser1][ag];
395 md->cU2[i] = groups->grpnr[egcUser2][ag];
400 if (groups->grpnr[egcQMMM] == nullptr ||
401 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
411 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
415 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
416 * For free-energy runs, these should be updated using update_mdatoms().
418 md->tmass = md->tmassA;
422 void update_mdatoms(t_mdatoms *md, real lambda)
424 if (md->nMassPerturbed && lambda != md->lambda)
426 real L1 = 1 - lambda;
428 /* Update masses of perturbed atoms for the change in lambda */
429 // cppcheck-suppress unreadVariable
430 int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
431 #pragma omp parallel for num_threads(nthreads) schedule(static)
432 for (int i = 0; i < md->nr; i++)
434 if (md->bPerturbed[i])
436 md->massT[i] = L1*md->massA[i] + lambda*md->massB[i];
437 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
438 * and their invmass does not depend on lambda.
440 if (md->invmass[i] > 1.1*ALMOST_ZERO)
442 md->invmass[i] = 1.0/md->massT[i];
443 for (int d = 0; d < DIM; d++)
445 if (md->invMassPerDim[i][d] > 1.1*ALMOST_ZERO)
447 md->invMassPerDim[i][d] = md->invmass[i];
454 /* Update the system mass for the change in lambda */
455 md->tmass = L1*md->tmassA + lambda*md->tmassB;