Further improve getDDGridSetup
[gromacs.git] / src / gromacs / mdlib / mdatoms.cpp
blob45202fe9edeec380e936b196d0e5e66c93991c47
1 /*
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,2018,2019, 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.
37 #include "gmxpre.h"
39 #include "mdatoms.h"
41 #include <cmath>
43 #include <memory>
45 #include "gromacs/ewald/pme.h"
46 #include "gromacs/gpu_utils/hostallocator.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/mdlib/gmx_omp_nthreads.h"
49 #include "gromacs/mdlib/qmmm.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/topology/mtop_lookup.h"
53 #include "gromacs/topology/mtop_util.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/smalloc.h"
58 #define ALMOST_ZERO 1e-30
60 namespace gmx
63 MDAtoms::MDAtoms()
64 : mdatoms_(nullptr)
68 MDAtoms::~MDAtoms()
70 if (mdatoms_ == nullptr)
72 return;
74 sfree(mdatoms_->massA);
75 sfree(mdatoms_->massB);
76 sfree(mdatoms_->massT);
77 gmx::AlignedAllocationPolicy::free(mdatoms_->invmass);
78 sfree(mdatoms_->invMassPerDim);
79 sfree(mdatoms_->typeA);
80 sfree(mdatoms_->chargeB);
81 sfree(mdatoms_->typeB);
82 sfree(mdatoms_->sqrt_c6A);
83 sfree(mdatoms_->sigmaA);
84 sfree(mdatoms_->sigma3A);
85 sfree(mdatoms_->sqrt_c6B);
86 sfree(mdatoms_->sigmaB);
87 sfree(mdatoms_->sigma3B);
88 sfree(mdatoms_->ptype);
89 sfree(mdatoms_->cTC);
90 sfree(mdatoms_->cENER);
91 sfree(mdatoms_->cACC);
92 sfree(mdatoms_->cFREEZE);
93 sfree(mdatoms_->cVCM);
94 sfree(mdatoms_->cORF);
95 sfree(mdatoms_->bPerturbed);
96 sfree(mdatoms_->cU1);
97 sfree(mdatoms_->cU2);
98 sfree(mdatoms_->bQM);
101 void MDAtoms::resize(int newSize)
103 chargeA_.resizeWithPadding(newSize);
104 mdatoms_->chargeA = chargeA_.data();
107 void MDAtoms::reserve(int newCapacity)
109 chargeA_.reserveWithPadding(newCapacity);
110 mdatoms_->chargeA = chargeA_.data();
113 std::unique_ptr<MDAtoms>
114 makeMDAtoms(FILE *fp, const gmx_mtop_t &mtop, const t_inputrec &ir,
115 const bool rankHasPmeGpuTask)
117 auto mdAtoms = std::make_unique<MDAtoms>();
118 // GPU transfers may want to use a suitable pinning mode.
119 if (rankHasPmeGpuTask)
121 changePinningPolicy(&mdAtoms->chargeA_, pme_get_pinning_policy());
123 t_mdatoms *md;
124 snew(md, 1);
125 mdAtoms->mdatoms_.reset(md);
127 md->nenergrp = mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size();
128 md->bVCMgrps = FALSE;
129 for (int i = 0; i < mtop.natoms; i++)
131 if (getGroupType(mtop.groups, SimulationAtomGroupType::MassCenterVelocityRemoval, i) > 0)
133 md->bVCMgrps = TRUE;
137 /* Determine the total system mass and perturbed atom counts */
138 double totalMassA = 0.0;
139 double totalMassB = 0.0;
141 md->haveVsites = FALSE;
142 gmx_mtop_atomloop_block_t aloop = gmx_mtop_atomloop_block_init(&mtop);
143 const t_atom *atom;
144 int nmol;
145 while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
147 totalMassA += nmol*atom->m;
148 totalMassB += nmol*atom->mB;
150 if (atom->ptype == eptVSite)
152 md->haveVsites = TRUE;
155 if (ir.efep != efepNO && PERTURBED(*atom))
157 md->nPerturbed++;
158 if (atom->mB != atom->m)
160 md->nMassPerturbed += nmol;
162 if (atom->qB != atom->q)
164 md->nChargePerturbed += nmol;
166 if (atom->typeB != atom->type)
168 md->nTypePerturbed += nmol;
173 md->tmassA = totalMassA;
174 md->tmassB = totalMassB;
176 if (ir.efep != efepNO && fp)
178 fprintf(fp,
179 "There are %d atoms and %d charges for free energy perturbation\n",
180 md->nPerturbed, md->nChargePerturbed);
183 md->havePartiallyFrozenAtoms = FALSE;
184 for (int g = 0; g < ir.opts.ngfrz; g++)
186 for (int d = YY; d < DIM; d++)
188 if (ir.opts.nFreeze[g][d] != ir.opts.nFreeze[g][XX])
190 md->havePartiallyFrozenAtoms = TRUE;
195 md->bOrires = (gmx_mtop_ftype_count(&mtop, F_ORIRES) != 0);
197 return mdAtoms;
200 } // namespace gmx
202 void atoms2md(const gmx_mtop_t *mtop, const t_inputrec *ir,
203 int nindex, const int *index,
204 int homenr,
205 gmx::MDAtoms *mdAtoms)
207 gmx_bool bLJPME;
208 const t_grpopts *opts;
209 int nthreads gmx_unused;
211 bLJPME = EVDW_PME(ir->vdwtype);
213 opts = &ir->opts;
215 const SimulationGroups &groups = mtop->groups;
217 auto md = mdAtoms->mdatoms();
218 /* nindex>=0 indicates DD where we use an index */
219 if (nindex >= 0)
221 md->nr = nindex;
223 else
225 md->nr = mtop->natoms;
228 if (md->nr > md->nalloc)
230 md->nalloc = over_alloc_dd(md->nr);
232 if (md->nMassPerturbed)
234 srenew(md->massA, md->nalloc);
235 srenew(md->massB, md->nalloc);
237 srenew(md->massT, md->nalloc);
238 /* The SIMD version of the integrator needs this aligned and padded.
239 * The padding needs to be with zeros, which we set later below.
241 gmx::AlignedAllocationPolicy::free(md->invmass);
242 md->invmass = new(gmx::AlignedAllocationPolicy::malloc((md->nalloc + GMX_REAL_MAX_SIMD_WIDTH)*sizeof(*md->invmass)))real;
243 srenew(md->invMassPerDim, md->nalloc);
244 // TODO eventually we will have vectors and just resize
245 // everything, but for now the semantics of md->nalloc being
246 // the capacity are preserved by keeping vectors within
247 // mdAtoms having the same properties as the other arrays.
248 mdAtoms->reserve(md->nalloc);
249 mdAtoms->resize(md->nr);
250 srenew(md->typeA, md->nalloc);
251 if (md->nPerturbed)
253 srenew(md->chargeB, md->nalloc);
254 srenew(md->typeB, md->nalloc);
256 if (bLJPME)
258 srenew(md->sqrt_c6A, md->nalloc);
259 srenew(md->sigmaA, md->nalloc);
260 srenew(md->sigma3A, md->nalloc);
261 if (md->nPerturbed)
263 srenew(md->sqrt_c6B, md->nalloc);
264 srenew(md->sigmaB, md->nalloc);
265 srenew(md->sigma3B, md->nalloc);
268 srenew(md->ptype, md->nalloc);
269 if (opts->ngtc > 1)
271 srenew(md->cTC, md->nalloc);
272 /* We always copy cTC with domain decomposition */
274 srenew(md->cENER, md->nalloc);
275 if (opts->ngacc > 1)
277 srenew(md->cACC, md->nalloc);
279 if (opts->nFreeze &&
280 (opts->ngfrz > 1 ||
281 opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
283 srenew(md->cFREEZE, md->nalloc);
285 if (md->bVCMgrps)
287 srenew(md->cVCM, md->nalloc);
289 if (md->bOrires)
291 srenew(md->cORF, md->nalloc);
293 if (md->nPerturbed)
295 srenew(md->bPerturbed, md->nalloc);
298 /* Note that these user t_mdatoms array pointers are NULL
299 * when there is only one group present.
300 * Therefore, when adding code, the user should use something like:
301 * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
303 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User1].empty())
305 srenew(md->cU1, md->nalloc);
307 if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User2].empty())
309 srenew(md->cU2, md->nalloc);
312 if (ir->bQMMM)
314 srenew(md->bQM, md->nalloc);
318 int molb = 0;
320 nthreads = gmx_omp_nthreads_get(emntDefault);
321 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
322 for (int i = 0; i < md->nr; i++)
326 int g, ag;
327 real mA, mB, fac;
328 real c6, c12;
330 if (index == nullptr)
332 ag = i;
334 else
336 ag = index[i];
338 const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
340 if (md->cFREEZE)
342 md->cFREEZE[i] = getGroupType(groups, SimulationAtomGroupType::Freeze, ag);
344 if (EI_ENERGY_MINIMIZATION(ir->eI))
346 /* Displacement is proportional to F, masses used for constraints */
347 mA = 1.0;
348 mB = 1.0;
350 else if (ir->eI == eiBD)
352 /* With BD the physical masses are irrelevant.
353 * To keep the code simple we use most of the normal MD code path
354 * for BD. Thus for constraining the masses should be proportional
355 * to the friction coefficient. We set the absolute value such that
356 * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
357 * Then if we set the (meaningless) velocity to v=dx/dt, we get the
358 * correct kinetic energy and temperature using the usual code path.
359 * Thus with BD v*dt will give the displacement and the reported
360 * temperature can signal bad integration (too large time step).
362 if (ir->bd_fric > 0)
364 mA = 0.5*ir->bd_fric*ir->delta_t;
365 mB = 0.5*ir->bd_fric*ir->delta_t;
367 else
369 /* The friction coefficient is mass/tau_t */
370 fac = ir->delta_t/opts->tau_t[md->cTC ?
371 groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag] : 0];
372 mA = 0.5*atom.m*fac;
373 mB = 0.5*atom.mB*fac;
376 else
378 mA = atom.m;
379 mB = atom.mB;
381 if (md->nMassPerturbed)
383 md->massA[i] = mA;
384 md->massB[i] = mB;
386 md->massT[i] = mA;
388 if (mA == 0.0)
390 md->invmass[i] = 0;
391 md->invMassPerDim[i][XX] = 0;
392 md->invMassPerDim[i][YY] = 0;
393 md->invMassPerDim[i][ZZ] = 0;
395 else if (md->cFREEZE)
397 g = md->cFREEZE[i];
398 GMX_ASSERT(opts->nFreeze != nullptr, "Must have freeze groups to initialize masses");
399 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
401 /* Set the mass of completely frozen particles to ALMOST_ZERO
402 * iso 0 to avoid div by zero in lincs or shake.
404 md->invmass[i] = ALMOST_ZERO;
406 else
408 /* Note: Partially frozen particles use the normal invmass.
409 * If such particles are constrained, the frozen dimensions
410 * should not be updated with the constrained coordinates.
412 md->invmass[i] = 1.0/mA;
414 for (int d = 0; d < DIM; d++)
416 md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0/mA);
419 else
421 md->invmass[i] = 1.0/mA;
422 for (int d = 0; d < DIM; d++)
424 md->invMassPerDim[i][d] = 1.0/mA;
428 md->chargeA[i] = atom.q;
429 md->typeA[i] = atom.type;
430 if (bLJPME)
432 c6 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c6;
433 c12 = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c12;
434 md->sqrt_c6A[i] = sqrt(c6);
435 if (c6 == 0.0 || c12 == 0)
437 md->sigmaA[i] = 1.0;
439 else
441 md->sigmaA[i] = gmx::sixthroot(c12/c6);
443 md->sigma3A[i] = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
445 if (md->nPerturbed)
447 md->bPerturbed[i] = PERTURBED(atom);
448 md->chargeB[i] = atom.qB;
449 md->typeB[i] = atom.typeB;
450 if (bLJPME)
452 c6 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c6;
453 c12 = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c12;
454 md->sqrt_c6B[i] = sqrt(c6);
455 if (c6 == 0.0 || c12 == 0)
457 md->sigmaB[i] = 1.0;
459 else
461 md->sigmaB[i] = gmx::sixthroot(c12/c6);
463 md->sigma3B[i] = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
466 md->ptype[i] = atom.ptype;
467 if (md->cTC)
469 md->cTC[i] = groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag];
471 md->cENER[i] = getGroupType(groups, SimulationAtomGroupType::EnergyOutput, ag);
472 if (md->cACC)
474 md->cACC[i] = groups.groupNumbers[SimulationAtomGroupType::Acceleration][ag];
476 if (md->cVCM)
478 md->cVCM[i] = groups.groupNumbers[SimulationAtomGroupType::MassCenterVelocityRemoval][ag];
480 if (md->cORF)
482 md->cORF[i] = getGroupType(groups, SimulationAtomGroupType::OrientationRestraintsFit, ag);
485 if (md->cU1)
487 md->cU1[i] = groups.groupNumbers[SimulationAtomGroupType::User1][ag];
489 if (md->cU2)
491 md->cU2[i] = groups.groupNumbers[SimulationAtomGroupType::User2][ag];
494 if (ir->bQMMM)
496 if (groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty() ||
497 groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][ag] < groups.groups[SimulationAtomGroupType::QuantumMechanics].size()-1)
499 md->bQM[i] = TRUE;
501 else
503 md->bQM[i] = FALSE;
507 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
510 if (md->nr > 0)
512 /* Pad invmass with 0 so a SIMD MD update does not change v and x */
513 for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
515 md->invmass[i] = 0;
519 md->homenr = homenr;
520 /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
521 * For free-energy runs, these should be updated using update_mdatoms().
523 md->tmass = md->tmassA;
524 md->lambda = 0;
527 void update_mdatoms(t_mdatoms *md, real lambda)
529 if (md->nMassPerturbed && lambda != md->lambda)
531 real L1 = 1 - lambda;
533 /* Update masses of perturbed atoms for the change in lambda */
534 int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
535 #pragma omp parallel for num_threads(nthreads) schedule(static)
536 for (int i = 0; i < md->nr; i++)
538 if (md->bPerturbed[i])
540 md->massT[i] = L1*md->massA[i] + lambda*md->massB[i];
541 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
542 * and their invmass does not depend on lambda.
544 if (md->invmass[i] > 1.1*ALMOST_ZERO)
546 md->invmass[i] = 1.0/md->massT[i];
547 for (int d = 0; d < DIM; d++)
549 if (md->invMassPerDim[i][d] > 1.1*ALMOST_ZERO)
551 md->invMassPerDim[i][d] = md->invmass[i];
558 /* Update the system mass for the change in lambda */
559 md->tmass = L1*md->tmassA + lambda*md->tmassB;
562 md->lambda = lambda;