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.
45 #include "gromacs/compat/make_unique.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/mdlib/gmx_omp_nthreads.h"
48 #include "gromacs/mdlib/qmmm.h"
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/topology/mtop_lookup.h"
52 #include "gromacs/topology/mtop_util.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/alignedallocator.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/smalloc.h"
58 #define ALMOST_ZERO 1e-30
63 MDAtoms::MDAtoms(HostAllocationPolicy policy)
64 : mdatoms_(nullptr), chargeA_(policy)
68 void MDAtoms::resize(int newSize)
70 chargeA_.resize(newSize);
71 mdatoms_->chargeA = chargeA_.data();
74 void MDAtoms::reserve(int newCapacity)
76 chargeA_.reserve(newCapacity);
77 mdatoms_->chargeA = chargeA_.data();
80 std::unique_ptr<MDAtoms>
81 makeMDAtoms(FILE *fp, const gmx_mtop_t &mtop, const t_inputrec &ir,
84 auto policy = (useGpuForPme ?
85 makeHostAllocationPolicyForGpu() :
86 HostAllocationPolicy());
87 auto mdAtoms = compat::make_unique<MDAtoms>(policy);
90 mdAtoms->mdatoms_.reset(md);
92 md->nenergrp = mtop.groups.grps[egcENER].nr;
93 md->bVCMgrps = (mtop.groups.grps[egcVCM].nr > 1);
95 /* Determine the total system mass and perturbed atom counts */
96 double totalMassA = 0.0;
97 double totalMassB = 0.0;
99 gmx_mtop_atomloop_block_t aloop = gmx_mtop_atomloop_block_init(&mtop);
102 while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
104 totalMassA += nmol*atom->m;
105 totalMassB += nmol*atom->mB;
107 if (ir.efep != efepNO && PERTURBED(*atom))
110 if (atom->mB != atom->m)
112 md->nMassPerturbed += nmol;
114 if (atom->qB != atom->q)
116 md->nChargePerturbed += nmol;
118 if (atom->typeB != atom->type)
120 md->nTypePerturbed += nmol;
125 md->tmassA = totalMassA;
126 md->tmassB = totalMassB;
128 if (ir.efep != efepNO && fp)
131 "There are %d atoms and %d charges for free energy perturbation\n",
132 md->nPerturbed, md->nChargePerturbed);
135 md->havePartiallyFrozenAtoms = FALSE;
136 for (int g = 0; g < ir.opts.ngfrz; g++)
138 for (int d = YY; d < DIM; d++)
140 if (ir.opts.nFreeze[d] != ir.opts.nFreeze[XX])
142 md->havePartiallyFrozenAtoms = TRUE;
147 md->bOrires = gmx_mtop_ftype_count(&mtop, F_ORIRES);
154 void atoms2md(const gmx_mtop_t *mtop, const t_inputrec *ir,
155 int nindex, const int *index,
157 gmx::MDAtoms *mdAtoms)
160 const t_grpopts *opts;
161 const gmx_groups_t *groups;
162 int nthreads gmx_unused;
164 bLJPME = EVDW_PME(ir->vdwtype);
168 groups = &mtop->groups;
170 auto md = mdAtoms->mdatoms();
171 /* nindex>=0 indicates DD where we use an index */
178 md->nr = mtop->natoms;
181 if (md->nr > md->nalloc)
183 md->nalloc = over_alloc_dd(md->nr);
185 if (md->nMassPerturbed)
187 srenew(md->massA, md->nalloc);
188 srenew(md->massB, md->nalloc);
190 srenew(md->massT, md->nalloc);
191 /* The SIMD version of the integrator needs this aligned and padded.
192 * The padding needs to be with zeros, which we set later below.
194 gmx::AlignedAllocationPolicy::free(md->invmass);
195 md->invmass = new(gmx::AlignedAllocationPolicy::malloc((md->nalloc + GMX_REAL_MAX_SIMD_WIDTH)*sizeof(*md->invmass)))real;
196 srenew(md->invMassPerDim, md->nalloc);
197 // TODO eventually we will have vectors and just resize
198 // everything, but for now the semantics of md->nalloc being
199 // the capacity are preserved by keeping vectors within
200 // mdAtoms having the same properties as the other arrays.
201 mdAtoms->reserve(md->nalloc);
202 mdAtoms->resize(md->nr);
203 srenew(md->typeA, md->nalloc);
206 srenew(md->chargeB, md->nalloc);
207 srenew(md->typeB, md->nalloc);
211 srenew(md->sqrt_c6A, md->nalloc);
212 srenew(md->sigmaA, md->nalloc);
213 srenew(md->sigma3A, md->nalloc);
216 srenew(md->sqrt_c6B, md->nalloc);
217 srenew(md->sigmaB, md->nalloc);
218 srenew(md->sigma3B, md->nalloc);
221 srenew(md->ptype, md->nalloc);
224 srenew(md->cTC, md->nalloc);
225 /* We always copy cTC with domain decomposition */
227 srenew(md->cENER, md->nalloc);
230 srenew(md->cACC, md->nalloc);
234 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
236 srenew(md->cFREEZE, md->nalloc);
240 srenew(md->cVCM, md->nalloc);
244 srenew(md->cORF, md->nalloc);
248 srenew(md->bPerturbed, md->nalloc);
251 /* Note that these user t_mdatoms array pointers are NULL
252 * when there is only one group present.
253 * Therefore, when adding code, the user should use something like:
254 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
256 if (mtop->groups.grpnr[egcUser1] != nullptr)
258 srenew(md->cU1, md->nalloc);
260 if (mtop->groups.grpnr[egcUser2] != nullptr)
262 srenew(md->cU2, md->nalloc);
267 srenew(md->bQM, md->nalloc);
273 // cppcheck-suppress unreadVariable
274 nthreads = gmx_omp_nthreads_get(emntDefault);
275 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
276 for (int i = 0; i < md->nr; i++)
284 if (index == nullptr)
292 const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
296 md->cFREEZE[i] = ggrpnr(groups, egcFREEZE, ag);
298 if (EI_ENERGY_MINIMIZATION(ir->eI))
300 /* Displacement is proportional to F, masses used for constraints */
304 else if (ir->eI == eiBD)
306 /* With BD the physical masses are irrelevant.
307 * To keep the code simple we use most of the normal MD code path
308 * for BD. Thus for constraining the masses should be proportional
309 * to the friction coefficient. We set the absolute value such that
310 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
311 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
312 * correct kinetic energy and temperature using the usual code path.
313 * Thus with BD v*dt will give the displacement and the reported
314 * temperature can signal bad integration (too large time step).
318 mA = 0.5*ir->bd_fric*ir->delta_t;
319 mB = 0.5*ir->bd_fric*ir->delta_t;
323 /* The friction coefficient is mass/tau_t */
324 fac = ir->delta_t/opts->tau_t[md->cTC ? groups->grpnr[egcTC][ag] : 0];
326 mB = 0.5*atom.mB*fac;
334 if (md->nMassPerturbed)
344 md->invMassPerDim[i][XX] = 0;
345 md->invMassPerDim[i][YY] = 0;
346 md->invMassPerDim[i][ZZ] = 0;
348 else if (md->cFREEZE)
351 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
353 /* Set the mass of completely frozen particles to ALMOST_ZERO
354 * iso 0 to avoid div by zero in lincs or shake.
356 md->invmass[i] = ALMOST_ZERO;
360 /* Note: Partially frozen particles use the normal invmass.
361 * If such particles are constrained, the frozen dimensions
362 * should not be updated with the constrained coordinates.
364 md->invmass[i] = 1.0/mA;
366 for (int d = 0; d < DIM; d++)
368 md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0/mA);
373 md->invmass[i] = 1.0/mA;
374 for (int d = 0; d < DIM; d++)
376 md->invMassPerDim[i][d] = 1.0/mA;
380 md->chargeA[i] = atom.q;
381 md->typeA[i] = atom.type;
384 c6 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c6;
385 c12 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c12;
386 md->sqrt_c6A[i] = sqrt(c6);
387 if (c6 == 0.0 || c12 == 0)
393 md->sigmaA[i] = gmx::sixthroot(c12/c6);
395 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
399 md->bPerturbed[i] = PERTURBED(atom);
400 md->chargeB[i] = atom.qB;
401 md->typeB[i] = atom.typeB;
404 c6 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c6;
405 c12 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c12;
406 md->sqrt_c6B[i] = sqrt(c6);
407 if (c6 == 0.0 || c12 == 0)
413 md->sigmaB[i] = gmx::sixthroot(c12/c6);
415 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
418 md->ptype[i] = atom.ptype;
421 md->cTC[i] = groups->grpnr[egcTC][ag];
423 md->cENER[i] = ggrpnr(groups, egcENER, ag);
426 md->cACC[i] = groups->grpnr[egcACC][ag];
430 md->cVCM[i] = groups->grpnr[egcVCM][ag];
434 md->cORF[i] = ggrpnr(groups, egcORFIT, ag);
439 md->cU1[i] = groups->grpnr[egcUser1][ag];
443 md->cU2[i] = groups->grpnr[egcUser2][ag];
448 if (groups->grpnr[egcQMMM] == nullptr ||
449 groups->grpnr[egcQMMM][ag] < groups->grps[egcQMMM].nr-1)
459 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
464 /* Pad invmass with 0 so a SIMD MD update does not change v and x */
465 for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
472 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
473 * For free-energy runs, these should be updated using update_mdatoms().
475 md->tmass = md->tmassA;
479 void update_mdatoms(t_mdatoms *md, real lambda)
481 if (md->nMassPerturbed && lambda != md->lambda)
483 real L1 = 1 - lambda;
485 /* Update masses of perturbed atoms for the change in lambda */
486 // cppcheck-suppress unreadVariable
487 int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
488 #pragma omp parallel for num_threads(nthreads) schedule(static)
489 for (int i = 0; i < md->nr; i++)
491 if (md->bPerturbed[i])
493 md->massT[i] = L1*md->massA[i] + lambda*md->massB[i];
494 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
495 * and their invmass does not depend on lambda.
497 if (md->invmass[i] > 1.1*ALMOST_ZERO)
499 md->invmass[i] = 1.0/md->massT[i];
500 for (int d = 0; d < DIM; d++)
502 if (md->invMassPerDim[i][d] > 1.1*ALMOST_ZERO)
504 md->invMassPerDim[i][d] = md->invmass[i];
511 /* Update the system mass for the change in lambda */
512 md->tmass = L1*md->tmassA + lambda*md->tmassB;