Reimplement constant acceleration groups
[alexxy/gromacs.git] / src / gromacs / mdlib / mdatoms.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "mdatoms.h"
41
42 #include <cmath>
43
44 #include <memory>
45
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"
58
59 #define ALMOST_ZERO 1e-30
60
61 namespace gmx
62 {
63
64 MDAtoms::MDAtoms() : mdatoms_(nullptr) {}
65
66 MDAtoms::~MDAtoms()
67 {
68     if (mdatoms_ == nullptr)
69     {
70         return;
71     }
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_->typeB);
79     /* mdatoms->chargeA and mdatoms->chargeB point at chargeA_.data()
80      * and chargeB_.data() respectively. They get freed automatically. */
81     sfree(mdatoms_->sqrt_c6A);
82     sfree(mdatoms_->sigmaA);
83     sfree(mdatoms_->sigma3A);
84     sfree(mdatoms_->sqrt_c6B);
85     sfree(mdatoms_->sigmaB);
86     sfree(mdatoms_->sigma3B);
87     sfree(mdatoms_->ptype);
88     sfree(mdatoms_->cTC);
89     sfree(mdatoms_->cENER);
90     sfree(mdatoms_->cACC);
91     sfree(mdatoms_->cFREEZE);
92     sfree(mdatoms_->cVCM);
93     sfree(mdatoms_->cORF);
94     sfree(mdatoms_->bPerturbed);
95     sfree(mdatoms_->cU1);
96     sfree(mdatoms_->cU2);
97 }
98
99 void MDAtoms::resizeChargeA(const int newSize)
100 {
101     chargeA_.resizeWithPadding(newSize);
102     mdatoms_->chargeA = chargeA_.data();
103 }
104
105 void MDAtoms::resizeChargeB(const int newSize)
106 {
107     chargeB_.resizeWithPadding(newSize);
108     mdatoms_->chargeB = chargeB_.data();
109 }
110
111 void MDAtoms::reserveChargeA(const int newCapacity)
112 {
113     chargeA_.reserveWithPadding(newCapacity);
114     mdatoms_->chargeA = chargeA_.data();
115 }
116
117 void MDAtoms::reserveChargeB(const int newCapacity)
118 {
119     chargeB_.reserveWithPadding(newCapacity);
120     mdatoms_->chargeB = chargeB_.data();
121 }
122
123 std::unique_ptr<MDAtoms> makeMDAtoms(FILE* fp, const gmx_mtop_t& mtop, const t_inputrec& ir, const bool rankHasPmeGpuTask)
124 {
125     auto mdAtoms = std::make_unique<MDAtoms>();
126     // GPU transfers may want to use a suitable pinning mode.
127     if (rankHasPmeGpuTask)
128     {
129         changePinningPolicy(&mdAtoms->chargeA_, pme_get_pinning_policy());
130         changePinningPolicy(&mdAtoms->chargeB_, pme_get_pinning_policy());
131     }
132     t_mdatoms* md;
133     snew(md, 1);
134     mdAtoms->mdatoms_.reset(md);
135
136     md->nenergrp = mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].size();
137     md->bVCMgrps = FALSE;
138     for (int i = 0; i < mtop.natoms; i++)
139     {
140         if (getGroupType(mtop.groups, SimulationAtomGroupType::MassCenterVelocityRemoval, i) > 0)
141         {
142             md->bVCMgrps = TRUE;
143         }
144     }
145
146     /* Determine the total system mass and perturbed atom counts */
147     double totalMassA = 0.0;
148     double totalMassB = 0.0;
149
150     md->haveVsites                  = FALSE;
151     gmx_mtop_atomloop_block_t aloop = gmx_mtop_atomloop_block_init(mtop);
152     const t_atom*             atom;
153     int                       nmol;
154     while (gmx_mtop_atomloop_block_next(aloop, &atom, &nmol))
155     {
156         totalMassA += nmol * atom->m;
157         totalMassB += nmol * atom->mB;
158
159         if (atom->ptype == ParticleType::VSite)
160         {
161             md->haveVsites = TRUE;
162         }
163
164         if (ir.efep != FreeEnergyPerturbationType::No && PERTURBED(*atom))
165         {
166             md->nPerturbed++;
167             if (atom->mB != atom->m)
168             {
169                 md->nMassPerturbed += nmol;
170             }
171             if (atom->qB != atom->q)
172             {
173                 md->nChargePerturbed += nmol;
174             }
175             if (atom->typeB != atom->type)
176             {
177                 md->nTypePerturbed += nmol;
178             }
179         }
180     }
181
182     md->tmassA = totalMassA;
183     md->tmassB = totalMassB;
184
185     if (ir.efep != FreeEnergyPerturbationType::No && fp)
186     {
187         fprintf(fp,
188                 "There are %d atoms and %d charges for free energy perturbation\n",
189                 md->nPerturbed,
190                 md->nChargePerturbed);
191     }
192
193     md->havePartiallyFrozenAtoms = FALSE;
194     for (int g = 0; g < ir.opts.ngfrz; g++)
195     {
196         for (int d = YY; d < DIM; d++)
197         {
198             if (ir.opts.nFreeze[g][d] != ir.opts.nFreeze[g][XX])
199             {
200                 md->havePartiallyFrozenAtoms = TRUE;
201             }
202         }
203     }
204
205     md->bOrires = (gmx_mtop_ftype_count(mtop, F_ORIRES) != 0);
206
207     return mdAtoms;
208 }
209
210 } // namespace gmx
211
212 void atoms2md(const gmx_mtop_t&  mtop,
213               const t_inputrec&  inputrec,
214               int                nindex,
215               gmx::ArrayRef<int> index,
216               int                homenr,
217               gmx::MDAtoms*      mdAtoms)
218 {
219     gmx_bool         bLJPME;
220     const t_grpopts* opts;
221     int nthreads     gmx_unused;
222
223     bLJPME = EVDW_PME(inputrec.vdwtype);
224
225     opts = &inputrec.opts;
226
227     const SimulationGroups& groups = mtop.groups;
228
229     auto* md = mdAtoms->mdatoms();
230     /* nindex>=0 indicates DD where we use an index */
231     if (nindex >= 0)
232     {
233         md->nr = nindex;
234     }
235     else
236     {
237         md->nr = mtop.natoms;
238     }
239
240     if (md->nr > md->nalloc)
241     {
242         md->nalloc = over_alloc_dd(md->nr);
243
244         if (md->nMassPerturbed)
245         {
246             srenew(md->massA, md->nalloc);
247             srenew(md->massB, md->nalloc);
248         }
249         srenew(md->massT, md->nalloc);
250         /* The SIMD version of the integrator needs this aligned and padded.
251          * The padding needs to be with zeros, which we set later below.
252          */
253         gmx::AlignedAllocationPolicy::free(md->invmass);
254         md->invmass = new (gmx::AlignedAllocationPolicy::malloc(
255                 (md->nalloc + GMX_REAL_MAX_SIMD_WIDTH) * sizeof(*md->invmass))) real;
256         srenew(md->invMassPerDim, md->nalloc);
257         // TODO eventually we will have vectors and just resize
258         // everything, but for now the semantics of md->nalloc being
259         // the capacity are preserved by keeping vectors within
260         // mdAtoms having the same properties as the other arrays.
261         mdAtoms->reserveChargeA(md->nalloc);
262         mdAtoms->resizeChargeA(md->nr);
263         if (md->nPerturbed > 0)
264         {
265             mdAtoms->reserveChargeB(md->nalloc);
266             mdAtoms->resizeChargeB(md->nr);
267         }
268         srenew(md->typeA, md->nalloc);
269         if (md->nPerturbed)
270         {
271             srenew(md->typeB, md->nalloc);
272         }
273         if (bLJPME)
274         {
275             srenew(md->sqrt_c6A, md->nalloc);
276             srenew(md->sigmaA, md->nalloc);
277             srenew(md->sigma3A, md->nalloc);
278             if (md->nPerturbed)
279             {
280                 srenew(md->sqrt_c6B, md->nalloc);
281                 srenew(md->sigmaB, md->nalloc);
282                 srenew(md->sigma3B, md->nalloc);
283             }
284         }
285         srenew(md->ptype, md->nalloc);
286         if (opts->ngtc > 1)
287         {
288             srenew(md->cTC, md->nalloc);
289             /* We always copy cTC with domain decomposition */
290         }
291         srenew(md->cENER, md->nalloc);
292         if (inputrec.useConstantAcceleration)
293         {
294             srenew(md->cACC, md->nalloc);
295         }
296         if (inputrecFrozenAtoms(&inputrec))
297         {
298             srenew(md->cFREEZE, md->nalloc);
299         }
300         if (md->bVCMgrps)
301         {
302             srenew(md->cVCM, md->nalloc);
303         }
304         if (md->bOrires)
305         {
306             srenew(md->cORF, md->nalloc);
307         }
308         if (md->nPerturbed)
309         {
310             srenew(md->bPerturbed, md->nalloc);
311         }
312
313         /* Note that these user t_mdatoms array pointers are NULL
314          * when there is only one group present.
315          * Therefore, when adding code, the user should use something like:
316          * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
317          */
318         if (!mtop.groups.groupNumbers[SimulationAtomGroupType::User1].empty())
319         {
320             srenew(md->cU1, md->nalloc);
321         }
322         if (!mtop.groups.groupNumbers[SimulationAtomGroupType::User2].empty())
323         {
324             srenew(md->cU2, md->nalloc);
325         }
326     }
327
328     int molb = 0;
329
330     nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Default);
331 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
332     for (int i = 0; i < md->nr; i++)
333     {
334         try
335         {
336             int  g, ag;
337             real mA, mB, fac;
338             real c6, c12;
339
340             if (index.empty())
341             {
342                 ag = i;
343             }
344             else
345             {
346                 ag = index[i];
347             }
348             const t_atom& atom = mtopGetAtomParameters(mtop, ag, &molb);
349
350             if (md->cFREEZE)
351             {
352                 md->cFREEZE[i] = getGroupType(groups, SimulationAtomGroupType::Freeze, ag);
353             }
354             if (EI_ENERGY_MINIMIZATION(inputrec.eI))
355             {
356                 /* Displacement is proportional to F, masses used for constraints */
357                 mA = 1.0;
358                 mB = 1.0;
359             }
360             else if (inputrec.eI == IntegrationAlgorithm::BD)
361             {
362                 /* With BD the physical masses are irrelevant.
363                  * To keep the code simple we use most of the normal MD code path
364                  * for BD. Thus for constraining the masses should be proportional
365                  * to the friction coefficient. We set the absolute value such that
366                  * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
367                  * Then if we set the (meaningless) velocity to v=dx/dt, we get the
368                  * correct kinetic energy and temperature using the usual code path.
369                  * Thus with BD v*dt will give the displacement and the reported
370                  * temperature can signal bad integration (too large time step).
371                  */
372                 if (inputrec.bd_fric > 0)
373                 {
374                     mA = 0.5 * inputrec.bd_fric * inputrec.delta_t;
375                     mB = 0.5 * inputrec.bd_fric * inputrec.delta_t;
376                 }
377                 else
378                 {
379                     /* The friction coefficient is mass/tau_t */
380                     fac = inputrec.delta_t
381                           / opts->tau_t[md->cTC ? groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag] : 0];
382                     mA = 0.5 * atom.m * fac;
383                     mB = 0.5 * atom.mB * fac;
384                 }
385             }
386             else
387             {
388                 mA = atom.m;
389                 mB = atom.mB;
390             }
391             if (md->nMassPerturbed)
392             {
393                 md->massA[i] = mA;
394                 md->massB[i] = mB;
395             }
396             md->massT[i] = mA;
397
398             if (mA == 0.0)
399             {
400                 md->invmass[i]           = 0;
401                 md->invMassPerDim[i][XX] = 0;
402                 md->invMassPerDim[i][YY] = 0;
403                 md->invMassPerDim[i][ZZ] = 0;
404             }
405             else if (md->cFREEZE)
406             {
407                 g = md->cFREEZE[i];
408                 GMX_ASSERT(opts->nFreeze != nullptr, "Must have freeze groups to initialize masses");
409                 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
410                 {
411                     /* Set the mass of completely frozen particles to ALMOST_ZERO
412                      * iso 0 to avoid div by zero in lincs or shake.
413                      */
414                     md->invmass[i] = ALMOST_ZERO;
415                 }
416                 else
417                 {
418                     /* Note: Partially frozen particles use the normal invmass.
419                      * If such particles are constrained, the frozen dimensions
420                      * should not be updated with the constrained coordinates.
421                      */
422                     md->invmass[i] = 1.0 / mA;
423                 }
424                 for (int d = 0; d < DIM; d++)
425                 {
426                     md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0 / mA);
427                 }
428             }
429             else
430             {
431                 md->invmass[i] = 1.0 / mA;
432                 for (int d = 0; d < DIM; d++)
433                 {
434                     md->invMassPerDim[i][d] = 1.0 / mA;
435                 }
436             }
437
438             md->chargeA[i] = atom.q;
439             md->typeA[i]   = atom.type;
440             if (bLJPME)
441             {
442                 c6  = mtop.ffparams.iparams[atom.type * (mtop.ffparams.atnr + 1)].lj.c6;
443                 c12 = mtop.ffparams.iparams[atom.type * (mtop.ffparams.atnr + 1)].lj.c12;
444                 md->sqrt_c6A[i] = std::sqrt(c6);
445                 if (c6 == 0.0 || c12 == 0)
446                 {
447                     md->sigmaA[i] = 1.0;
448                 }
449                 else
450                 {
451                     md->sigmaA[i] = gmx::sixthroot(c12 / c6);
452                 }
453                 md->sigma3A[i] = 1 / (md->sigmaA[i] * md->sigmaA[i] * md->sigmaA[i]);
454             }
455             if (md->nPerturbed)
456             {
457                 md->bPerturbed[i] = PERTURBED(atom);
458                 md->chargeB[i]    = atom.qB;
459                 md->typeB[i]      = atom.typeB;
460                 if (bLJPME)
461                 {
462                     c6  = mtop.ffparams.iparams[atom.typeB * (mtop.ffparams.atnr + 1)].lj.c6;
463                     c12 = mtop.ffparams.iparams[atom.typeB * (mtop.ffparams.atnr + 1)].lj.c12;
464                     md->sqrt_c6B[i] = std::sqrt(c6);
465                     if (c6 == 0.0 || c12 == 0)
466                     {
467                         md->sigmaB[i] = 1.0;
468                     }
469                     else
470                     {
471                         md->sigmaB[i] = gmx::sixthroot(c12 / c6);
472                     }
473                     md->sigma3B[i] = 1 / (md->sigmaB[i] * md->sigmaB[i] * md->sigmaB[i]);
474                 }
475             }
476             md->ptype[i] = atom.ptype;
477             if (md->cTC)
478             {
479                 md->cTC[i] = groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag];
480             }
481             md->cENER[i] = getGroupType(groups, SimulationAtomGroupType::EnergyOutput, ag);
482             if (md->cACC)
483             {
484                 md->cACC[i] = groups.groupNumbers[SimulationAtomGroupType::Acceleration][ag];
485             }
486             if (md->cVCM)
487             {
488                 md->cVCM[i] = groups.groupNumbers[SimulationAtomGroupType::MassCenterVelocityRemoval][ag];
489             }
490             if (md->cORF)
491             {
492                 md->cORF[i] = getGroupType(groups, SimulationAtomGroupType::OrientationRestraintsFit, ag);
493             }
494
495             if (md->cU1)
496             {
497                 md->cU1[i] = groups.groupNumbers[SimulationAtomGroupType::User1][ag];
498             }
499             if (md->cU2)
500             {
501                 md->cU2[i] = groups.groupNumbers[SimulationAtomGroupType::User2][ag];
502             }
503         }
504         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
505     }
506
507     if (md->nr > 0)
508     {
509         /* Pad invmass with 0 so a SIMD MD update does not change v and x */
510         for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
511         {
512             md->invmass[i] = 0;
513         }
514     }
515
516     md->homenr = homenr;
517     /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
518      * For free-energy runs, these should be updated using update_mdatoms().
519      */
520     md->tmass  = md->tmassA;
521     md->lambda = 0;
522 }
523
524 void update_mdatoms(t_mdatoms* md, real lambda)
525 {
526     if (md->nMassPerturbed && lambda != md->lambda)
527     {
528         real L1 = 1 - lambda;
529
530         /* Update masses of perturbed atoms for the change in lambda */
531         int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Default);
532 #pragma omp parallel for num_threads(nthreads) schedule(static)
533         for (int i = 0; i < md->nr; i++)
534         {
535             if (md->bPerturbed[i])
536             {
537                 md->massT[i] = L1 * md->massA[i] + lambda * md->massB[i];
538                 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
539                  * and their invmass does not depend on lambda.
540                  */
541                 if (md->invmass[i] > 1.1 * ALMOST_ZERO)
542                 {
543                     md->invmass[i] = 1.0 / md->massT[i];
544                     for (int d = 0; d < DIM; d++)
545                     {
546                         if (md->invMassPerDim[i][d] > 1.1 * ALMOST_ZERO)
547                         {
548                             md->invMassPerDim[i][d] = md->invmass[i];
549                         }
550                     }
551                 }
552             }
553         }
554
555         /* Update the system mass for the change in lambda */
556         md->tmass = L1 * md->tmassA + lambda * md->tmassB;
557     }
558
559     md->lambda = lambda;
560 }