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 by the GROMACS development team.
7 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gromacs/ewald/pme.h"
47 #include "gromacs/gpu_utils/hostallocator.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/mdlib/gmx_omp_nthreads.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/topology/mtop_lookup.h"
54 #include "gromacs/topology/mtop_util.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/smalloc.h"
59 #define ALMOST_ZERO 1e-30
64 MDAtoms::MDAtoms() : mdatoms_(nullptr) {}
68 if (mdatoms_ == nullptr)
72 sfree(mdatoms_->massA);
73 sfree(mdatoms_->massB);
74 sfree(mdatoms_->massT);
75 gmx::AlignedAllocationPolicy::free(mdatoms_->invmass);
76 sfree(mdatoms_->invMassPerDim);
77 sfree(mdatoms_->typeA);
78 sfree(mdatoms_->chargeB);
79 sfree(mdatoms_->typeB);
80 sfree(mdatoms_->sqrt_c6A);
81 sfree(mdatoms_->sigmaA);
82 sfree(mdatoms_->sigma3A);
83 sfree(mdatoms_->sqrt_c6B);
84 sfree(mdatoms_->sigmaB);
85 sfree(mdatoms_->sigma3B);
86 sfree(mdatoms_->ptype);
88 sfree(mdatoms_->cENER);
89 sfree(mdatoms_->cACC);
90 sfree(mdatoms_->cFREEZE);
91 sfree(mdatoms_->cVCM);
92 sfree(mdatoms_->cORF);
93 sfree(mdatoms_->bPerturbed);
98 void MDAtoms::resize(int newSize)
100 chargeA_.resizeWithPadding(newSize);
101 mdatoms_->chargeA = chargeA_.data();
104 void MDAtoms::reserve(int newCapacity)
106 chargeA_.reserveWithPadding(newCapacity);
107 mdatoms_->chargeA = chargeA_.data();
110 std::unique_ptr<MDAtoms> makeMDAtoms(FILE* fp, const gmx_mtop_t& mtop, const t_inputrec& ir, const bool rankHasPmeGpuTask)
112 auto mdAtoms = std::make_unique<MDAtoms>();
113 // GPU transfers may want to use a suitable pinning mode.
114 if (rankHasPmeGpuTask)
116 changePinningPolicy(&mdAtoms->chargeA_, pme_get_pinning_policy());
120 mdAtoms->mdatoms_.reset(md);
122 md->nenergrp = mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size();
123 md->bVCMgrps = FALSE;
124 for (int i = 0; i < mtop.natoms; i++)
126 if (getGroupType(mtop.groups, SimulationAtomGroupType::MassCenterVelocityRemoval, i) > 0)
132 /* Determine the total system mass and perturbed atom counts */
133 double totalMassA = 0.0;
134 double totalMassB = 0.0;
136 md->haveVsites = FALSE;
137 gmx_mtop_atomloop_block_t aloop = gmx_mtop_atomloop_block_init(&mtop);
140 while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
142 totalMassA += nmol * atom->m;
143 totalMassB += nmol * atom->mB;
145 if (atom->ptype == eptVSite)
147 md->haveVsites = TRUE;
150 if (ir.efep != efepNO && PERTURBED(*atom))
153 if (atom->mB != atom->m)
155 md->nMassPerturbed += nmol;
157 if (atom->qB != atom->q)
159 md->nChargePerturbed += nmol;
161 if (atom->typeB != atom->type)
163 md->nTypePerturbed += nmol;
168 md->tmassA = totalMassA;
169 md->tmassB = totalMassB;
171 if (ir.efep != efepNO && fp)
173 fprintf(fp, "There are %d atoms and %d charges for free energy perturbation\n",
174 md->nPerturbed, md->nChargePerturbed);
177 md->havePartiallyFrozenAtoms = FALSE;
178 for (int g = 0; g < ir.opts.ngfrz; g++)
180 for (int d = YY; d < DIM; d++)
182 if (ir.opts.nFreeze[g][d] != ir.opts.nFreeze[g][XX])
184 md->havePartiallyFrozenAtoms = TRUE;
189 md->bOrires = (gmx_mtop_ftype_count(&mtop, F_ORIRES) != 0);
196 void atoms2md(const gmx_mtop_t* mtop,
197 const t_inputrec* ir,
199 gmx::ArrayRef<int> index,
201 gmx::MDAtoms* mdAtoms)
204 const t_grpopts* opts;
205 int nthreads gmx_unused;
207 bLJPME = EVDW_PME(ir->vdwtype);
211 const SimulationGroups& groups = mtop->groups;
213 auto md = mdAtoms->mdatoms();
214 /* nindex>=0 indicates DD where we use an index */
221 md->nr = mtop->natoms;
224 if (md->nr > md->nalloc)
226 md->nalloc = over_alloc_dd(md->nr);
228 if (md->nMassPerturbed)
230 srenew(md->massA, md->nalloc);
231 srenew(md->massB, md->nalloc);
233 srenew(md->massT, md->nalloc);
234 /* The SIMD version of the integrator needs this aligned and padded.
235 * The padding needs to be with zeros, which we set later below.
237 gmx::AlignedAllocationPolicy::free(md->invmass);
238 md->invmass = new (gmx::AlignedAllocationPolicy::malloc(
239 (md->nalloc + GMX_REAL_MAX_SIMD_WIDTH) * sizeof(*md->invmass))) real;
240 srenew(md->invMassPerDim, md->nalloc);
241 // TODO eventually we will have vectors and just resize
242 // everything, but for now the semantics of md->nalloc being
243 // the capacity are preserved by keeping vectors within
244 // mdAtoms having the same properties as the other arrays.
245 mdAtoms->reserve(md->nalloc);
246 mdAtoms->resize(md->nr);
247 srenew(md->typeA, md->nalloc);
250 srenew(md->chargeB, md->nalloc);
251 srenew(md->typeB, md->nalloc);
255 srenew(md->sqrt_c6A, md->nalloc);
256 srenew(md->sigmaA, md->nalloc);
257 srenew(md->sigma3A, md->nalloc);
260 srenew(md->sqrt_c6B, md->nalloc);
261 srenew(md->sigmaB, md->nalloc);
262 srenew(md->sigma3B, md->nalloc);
265 srenew(md->ptype, md->nalloc);
268 srenew(md->cTC, md->nalloc);
269 /* We always copy cTC with domain decomposition */
271 srenew(md->cENER, md->nalloc);
274 srenew(md->cACC, md->nalloc);
277 && (opts->ngfrz > 1 || opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
279 srenew(md->cFREEZE, md->nalloc);
283 srenew(md->cVCM, md->nalloc);
287 srenew(md->cORF, md->nalloc);
291 srenew(md->bPerturbed, md->nalloc);
294 /* Note that these user t_mdatoms array pointers are NULL
295 * when there is only one group present.
296 * Therefore, when adding code, the user should use something like:
297 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
299 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User1].empty())
301 srenew(md->cU1, md->nalloc);
303 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User2].empty())
305 srenew(md->cU2, md->nalloc);
311 nthreads = gmx_omp_nthreads_get(emntDefault);
312 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
313 for (int i = 0; i < md->nr; i++)
329 const t_atom& atom = mtopGetAtomParameters(mtop, ag, &molb);
333 md->cFREEZE[i] = getGroupType(groups, SimulationAtomGroupType::Freeze, ag);
335 if (EI_ENERGY_MINIMIZATION(ir->eI))
337 /* Displacement is proportional to F, masses used for constraints */
341 else if (ir->eI == eiBD)
343 /* With BD the physical masses are irrelevant.
344 * To keep the code simple we use most of the normal MD code path
345 * for BD. Thus for constraining the masses should be proportional
346 * to the friction coefficient. We set the absolute value such that
347 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
348 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
349 * correct kinetic energy and temperature using the usual code path.
350 * Thus with BD v*dt will give the displacement and the reported
351 * temperature can signal bad integration (too large time step).
355 mA = 0.5 * ir->bd_fric * ir->delta_t;
356 mB = 0.5 * ir->bd_fric * ir->delta_t;
360 /* The friction coefficient is mass/tau_t */
362 / opts->tau_t[md->cTC ? groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag] : 0];
363 mA = 0.5 * atom.m * fac;
364 mB = 0.5 * atom.mB * fac;
372 if (md->nMassPerturbed)
382 md->invMassPerDim[i][XX] = 0;
383 md->invMassPerDim[i][YY] = 0;
384 md->invMassPerDim[i][ZZ] = 0;
386 else if (md->cFREEZE)
389 GMX_ASSERT(opts->nFreeze != nullptr,
390 "Must have freeze groups to initialize masses");
391 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
393 /* Set the mass of completely frozen particles to ALMOST_ZERO
394 * iso 0 to avoid div by zero in lincs or shake.
396 md->invmass[i] = ALMOST_ZERO;
400 /* Note: Partially frozen particles use the normal invmass.
401 * If such particles are constrained, the frozen dimensions
402 * should not be updated with the constrained coordinates.
404 md->invmass[i] = 1.0 / mA;
406 for (int d = 0; d < DIM; d++)
408 md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0 / mA);
413 md->invmass[i] = 1.0 / mA;
414 for (int d = 0; d < DIM; d++)
416 md->invMassPerDim[i][d] = 1.0 / mA;
420 md->chargeA[i] = atom.q;
421 md->typeA[i] = atom.type;
424 c6 = mtop->ffparams.iparams[atom.type * (mtop->ffparams.atnr + 1)].lj.c6;
425 c12 = mtop->ffparams.iparams[atom.type * (mtop->ffparams.atnr + 1)].lj.c12;
426 md->sqrt_c6A[i] = sqrt(c6);
427 if (c6 == 0.0 || c12 == 0)
433 md->sigmaA[i] = gmx::sixthroot(c12 / c6);
435 md->sigma3A[i] = 1 / (md->sigmaA[i] * md->sigmaA[i] * md->sigmaA[i]);
439 md->bPerturbed[i] = PERTURBED(atom);
440 md->chargeB[i] = atom.qB;
441 md->typeB[i] = atom.typeB;
444 c6 = mtop->ffparams.iparams[atom.typeB * (mtop->ffparams.atnr + 1)].lj.c6;
445 c12 = mtop->ffparams.iparams[atom.typeB * (mtop->ffparams.atnr + 1)].lj.c12;
446 md->sqrt_c6B[i] = sqrt(c6);
447 if (c6 == 0.0 || c12 == 0)
453 md->sigmaB[i] = gmx::sixthroot(c12 / c6);
455 md->sigma3B[i] = 1 / (md->sigmaB[i] * md->sigmaB[i] * md->sigmaB[i]);
458 md->ptype[i] = atom.ptype;
461 md->cTC[i] = groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag];
463 md->cENER[i] = getGroupType(groups, SimulationAtomGroupType::EnergyOutput, ag);
466 md->cACC[i] = groups.groupNumbers[SimulationAtomGroupType::Acceleration][ag];
470 md->cVCM[i] = groups.groupNumbers[SimulationAtomGroupType::MassCenterVelocityRemoval][ag];
474 md->cORF[i] = getGroupType(groups, SimulationAtomGroupType::OrientationRestraintsFit, ag);
479 md->cU1[i] = groups.groupNumbers[SimulationAtomGroupType::User1][ag];
483 md->cU2[i] = groups.groupNumbers[SimulationAtomGroupType::User2][ag];
486 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
491 /* Pad invmass with 0 so a SIMD MD update does not change v and x */
492 for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
499 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
500 * For free-energy runs, these should be updated using update_mdatoms().
502 md->tmass = md->tmassA;
506 void update_mdatoms(t_mdatoms* md, real lambda)
508 if (md->nMassPerturbed && lambda != md->lambda)
510 real L1 = 1 - lambda;
512 /* Update masses of perturbed atoms for the change in lambda */
513 int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
514 #pragma omp parallel for num_threads(nthreads) schedule(static)
515 for (int i = 0; i < md->nr; i++)
517 if (md->bPerturbed[i])
519 md->massT[i] = L1 * md->massA[i] + lambda * md->massB[i];
520 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
521 * and their invmass does not depend on lambda.
523 if (md->invmass[i] > 1.1 * ALMOST_ZERO)
525 md->invmass[i] = 1.0 / md->massT[i];
526 for (int d = 0; d < DIM; d++)
528 if (md->invMassPerDim[i][d] > 1.1 * ALMOST_ZERO)
530 md->invMassPerDim[i][d] = md->invmass[i];
537 /* Update the system mass for the change in lambda */
538 md->tmass = L1 * md->tmassA + lambda * md->tmassB;