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