84d73f4fd386c4afe3457cee02730aa0156b7953
[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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "mdatoms.h"
40
41 #include <cmath>
42
43 #include <memory>
44
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"
57
58 #define ALMOST_ZERO 1e-30
59
60 namespace gmx
61 {
62
63 MDAtoms::MDAtoms()
64     : mdatoms_(nullptr)
65 {
66 }
67
68 MDAtoms::~MDAtoms()
69 {
70     if (mdatoms_ == nullptr)
71     {
72         return;
73     }
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);
99 }
100
101 void MDAtoms::resize(int newSize)
102 {
103     chargeA_.resizeWithPadding(newSize);
104     mdatoms_->chargeA = chargeA_.data();
105 }
106
107 void MDAtoms::reserve(int newCapacity)
108 {
109     chargeA_.reserveWithPadding(newCapacity);
110     mdatoms_->chargeA = chargeA_.data();
111 }
112
113 std::unique_ptr<MDAtoms>
114 makeMDAtoms(FILE *fp, const gmx_mtop_t &mtop, const t_inputrec &ir,
115             const bool rankHasPmeGpuTask)
116 {
117     auto       mdAtoms = std::make_unique<MDAtoms>();
118     // GPU transfers may want to use a suitable pinning mode.
119     if (rankHasPmeGpuTask)
120     {
121         changePinningPolicy(&mdAtoms->chargeA_, pme_get_pinning_policy());
122     }
123     t_mdatoms *md;
124     snew(md, 1);
125     mdAtoms->mdatoms_.reset(md);
126
127     md->nenergrp = mtop.groups.groups[SimulationAtomGroupType::EnergyOutput].nr;
128     md->bVCMgrps = FALSE;
129     for (int i = 0; i < mtop.natoms; i++)
130     {
131         if (getGroupType(mtop.groups, SimulationAtomGroupType::MassCenterVelocityRemoval, i) > 0)
132         {
133             md->bVCMgrps = TRUE;
134         }
135     }
136
137     /* Determine the total system mass and perturbed atom counts */
138     double                     totalMassA = 0.0;
139     double                     totalMassB = 0.0;
140
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))
146     {
147         totalMassA += nmol*atom->m;
148         totalMassB += nmol*atom->mB;
149
150         if (atom->ptype == eptVSite)
151         {
152             md->haveVsites = TRUE;
153         }
154
155         if (ir.efep != efepNO && PERTURBED(*atom))
156         {
157             md->nPerturbed++;
158             if (atom->mB != atom->m)
159             {
160                 md->nMassPerturbed += nmol;
161             }
162             if (atom->qB != atom->q)
163             {
164                 md->nChargePerturbed += nmol;
165             }
166             if (atom->typeB != atom->type)
167             {
168                 md->nTypePerturbed += nmol;
169             }
170         }
171     }
172
173     md->tmassA = totalMassA;
174     md->tmassB = totalMassB;
175
176     if (ir.efep != efepNO && fp)
177     {
178         fprintf(fp,
179                 "There are %d atoms and %d charges for free energy perturbation\n",
180                 md->nPerturbed, md->nChargePerturbed);
181     }
182
183     md->havePartiallyFrozenAtoms = FALSE;
184     for (int g = 0; g < ir.opts.ngfrz; g++)
185     {
186         for (int d = YY; d < DIM; d++)
187         {
188             if (ir.opts.nFreeze[g][d] != ir.opts.nFreeze[g][XX])
189             {
190                 md->havePartiallyFrozenAtoms = TRUE;
191             }
192         }
193     }
194
195     md->bOrires = (gmx_mtop_ftype_count(&mtop, F_ORIRES) != 0);
196
197     return mdAtoms;
198 }
199
200 }  // namespace gmx
201
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)
206 {
207     gmx_bool              bLJPME;
208     const t_grpopts      *opts;
209     int                   nthreads gmx_unused;
210
211     bLJPME = EVDW_PME(ir->vdwtype);
212
213     opts = &ir->opts;
214
215     const SimulationGroups    &groups = mtop->groups;
216
217     auto                       md = mdAtoms->mdatoms();
218     /* nindex>=0 indicates DD where we use an index */
219     if (nindex >= 0)
220     {
221         md->nr = nindex;
222     }
223     else
224     {
225         md->nr = mtop->natoms;
226     }
227
228     if (md->nr > md->nalloc)
229     {
230         md->nalloc = over_alloc_dd(md->nr);
231
232         if (md->nMassPerturbed)
233         {
234             srenew(md->massA, md->nalloc);
235             srenew(md->massB, md->nalloc);
236         }
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.
240          */
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)
252         {
253             srenew(md->chargeB, md->nalloc);
254             srenew(md->typeB, md->nalloc);
255         }
256         if (bLJPME)
257         {
258             srenew(md->sqrt_c6A, md->nalloc);
259             srenew(md->sigmaA, md->nalloc);
260             srenew(md->sigma3A, md->nalloc);
261             if (md->nPerturbed)
262             {
263                 srenew(md->sqrt_c6B, md->nalloc);
264                 srenew(md->sigmaB, md->nalloc);
265                 srenew(md->sigma3B, md->nalloc);
266             }
267         }
268         srenew(md->ptype, md->nalloc);
269         if (opts->ngtc > 1)
270         {
271             srenew(md->cTC, md->nalloc);
272             /* We always copy cTC with domain decomposition */
273         }
274         srenew(md->cENER, md->nalloc);
275         if (opts->ngacc > 1)
276         {
277             srenew(md->cACC, md->nalloc);
278         }
279         if (opts->nFreeze &&
280             (opts->ngfrz > 1 ||
281              opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
282         {
283             srenew(md->cFREEZE, md->nalloc);
284         }
285         if (md->bVCMgrps)
286         {
287             srenew(md->cVCM, md->nalloc);
288         }
289         if (md->bOrires)
290         {
291             srenew(md->cORF, md->nalloc);
292         }
293         if (md->nPerturbed)
294         {
295             srenew(md->bPerturbed, md->nalloc);
296         }
297
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])
302          */
303         if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User1].empty())
304         {
305             srenew(md->cU1, md->nalloc);
306         }
307         if (!mtop->groups.groupNumbers[SimulationAtomGroupType::User2].empty())
308         {
309             srenew(md->cU2, md->nalloc);
310         }
311
312         if (ir->bQMMM)
313         {
314             srenew(md->bQM, md->nalloc);
315         }
316     }
317
318     int molb          = 0;
319
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++)
323     {
324         try
325         {
326             int      g, ag;
327             real     mA, mB, fac;
328             real     c6, c12;
329
330             if (index == nullptr)
331             {
332                 ag = i;
333             }
334             else
335             {
336                 ag = index[i];
337             }
338             const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
339
340             if (md->cFREEZE)
341             {
342                 md->cFREEZE[i] = getGroupType(groups, SimulationAtomGroupType::Freeze, ag);
343             }
344             if (EI_ENERGY_MINIMIZATION(ir->eI))
345             {
346                 /* Displacement is proportional to F, masses used for constraints */
347                 mA = 1.0;
348                 mB = 1.0;
349             }
350             else if (ir->eI == eiBD)
351             {
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).
361                  */
362                 if (ir->bd_fric > 0)
363                 {
364                     mA = 0.5*ir->bd_fric*ir->delta_t;
365                     mB = 0.5*ir->bd_fric*ir->delta_t;
366                 }
367                 else
368                 {
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;
374                 }
375             }
376             else
377             {
378                 mA = atom.m;
379                 mB = atom.mB;
380             }
381             if (md->nMassPerturbed)
382             {
383                 md->massA[i]  = mA;
384                 md->massB[i]  = mB;
385             }
386             md->massT[i]    = mA;
387
388             if (mA == 0.0)
389             {
390                 md->invmass[i]           = 0;
391                 md->invMassPerDim[i][XX] = 0;
392                 md->invMassPerDim[i][YY] = 0;
393                 md->invMassPerDim[i][ZZ] = 0;
394             }
395             else if (md->cFREEZE)
396             {
397                 g = md->cFREEZE[i];
398                 if (opts->nFreeze[g][XX] && opts->nFreeze[g][YY] && opts->nFreeze[g][ZZ])
399                 {
400                     /* Set the mass of completely frozen particles to ALMOST_ZERO
401                      * iso 0 to avoid div by zero in lincs or shake.
402                      */
403                     md->invmass[i]  = ALMOST_ZERO;
404                 }
405                 else
406                 {
407                     /* Note: Partially frozen particles use the normal invmass.
408                      * If such particles are constrained, the frozen dimensions
409                      * should not be updated with the constrained coordinates.
410                      */
411                     md->invmass[i]  = 1.0/mA;
412                 }
413                 for (int d = 0; d < DIM; d++)
414                 {
415                     md->invMassPerDim[i][d] = (opts->nFreeze[g][d] ? 0 : 1.0/mA);
416                 }
417             }
418             else
419             {
420                 md->invmass[i]  = 1.0/mA;
421                 for (int d = 0; d < DIM; d++)
422                 {
423                     md->invMassPerDim[i][d] = 1.0/mA;
424                 }
425             }
426
427             md->chargeA[i]      = atom.q;
428             md->typeA[i]        = atom.type;
429             if (bLJPME)
430             {
431                 c6                = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c6;
432                 c12               = mtop->ffparams.iparams[atom.type*(mtop->ffparams.atnr+1)].lj.c12;
433                 md->sqrt_c6A[i]   = sqrt(c6);
434                 if (c6 == 0.0 || c12 == 0)
435                 {
436                     md->sigmaA[i] = 1.0;
437                 }
438                 else
439                 {
440                     md->sigmaA[i] = gmx::sixthroot(c12/c6);
441                 }
442                 md->sigma3A[i]    = 1/(md->sigmaA[i]*md->sigmaA[i]*md->sigmaA[i]);
443             }
444             if (md->nPerturbed)
445             {
446                 md->bPerturbed[i] = PERTURBED(atom);
447                 md->chargeB[i]    = atom.qB;
448                 md->typeB[i]      = atom.typeB;
449                 if (bLJPME)
450                 {
451                     c6                = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c6;
452                     c12               = mtop->ffparams.iparams[atom.typeB*(mtop->ffparams.atnr+1)].lj.c12;
453                     md->sqrt_c6B[i]   = sqrt(c6);
454                     if (c6 == 0.0 || c12 == 0)
455                     {
456                         md->sigmaB[i] = 1.0;
457                     }
458                     else
459                     {
460                         md->sigmaB[i] = gmx::sixthroot(c12/c6);
461                     }
462                     md->sigma3B[i]    = 1/(md->sigmaB[i]*md->sigmaB[i]*md->sigmaB[i]);
463                 }
464             }
465             md->ptype[i]    = atom.ptype;
466             if (md->cTC)
467             {
468                 md->cTC[i]    = groups.groupNumbers[SimulationAtomGroupType::TemperatureCoupling][ag];
469             }
470             md->cENER[i]    = getGroupType(groups, SimulationAtomGroupType::EnergyOutput, ag);
471             if (md->cACC)
472             {
473                 md->cACC[i]   = groups.groupNumbers[SimulationAtomGroupType::Acceleration][ag];
474             }
475             if (md->cVCM)
476             {
477                 md->cVCM[i]       = groups.groupNumbers[SimulationAtomGroupType::MassCenterVelocityRemoval][ag];
478             }
479             if (md->cORF)
480             {
481                 md->cORF[i]       = getGroupType(groups, SimulationAtomGroupType::OrientationRestraintsFit, ag);
482             }
483
484             if (md->cU1)
485             {
486                 md->cU1[i]        = groups.groupNumbers[SimulationAtomGroupType::User1][ag];
487             }
488             if (md->cU2)
489             {
490                 md->cU2[i]        = groups.groupNumbers[SimulationAtomGroupType::User2][ag];
491             }
492
493             if (ir->bQMMM)
494             {
495                 if (groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics].empty() ||
496                     groups.groupNumbers[SimulationAtomGroupType::QuantumMechanics][ag] < groups.groups[SimulationAtomGroupType::QuantumMechanics].nr-1)
497                 {
498                     md->bQM[i]      = TRUE;
499                 }
500                 else
501                 {
502                     md->bQM[i]      = FALSE;
503                 }
504             }
505         }
506         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
507     }
508
509     if (md->nr > 0)
510     {
511         /* Pad invmass with 0 so a SIMD MD update does not change v and x */
512         for (int i = md->nr; i < md->nr + GMX_REAL_MAX_SIMD_WIDTH; i++)
513         {
514             md->invmass[i] = 0;
515         }
516     }
517
518     md->homenr = homenr;
519     /* We set mass, invmass, invMassPerDim and tmass for lambda=0.
520      * For free-energy runs, these should be updated using update_mdatoms().
521      */
522     md->tmass  = md->tmassA;
523     md->lambda = 0;
524 }
525
526 void update_mdatoms(t_mdatoms *md, real lambda)
527 {
528     if (md->nMassPerturbed && lambda != md->lambda)
529     {
530         real L1 = 1 - lambda;
531
532         /* Update masses of perturbed atoms for the change in lambda */
533         int gmx_unused nthreads = gmx_omp_nthreads_get(emntDefault);
534 #pragma omp parallel for num_threads(nthreads) schedule(static)
535         for (int i = 0; i < md->nr; i++)
536         {
537             if (md->bPerturbed[i])
538             {
539                 md->massT[i] = L1*md->massA[i] + lambda*md->massB[i];
540                 /* Atoms with invmass 0 or ALMOST_ZERO are massless or frozen
541                  * and their invmass does not depend on lambda.
542                  */
543                 if (md->invmass[i] > 1.1*ALMOST_ZERO)
544                 {
545                     md->invmass[i] = 1.0/md->massT[i];
546                     for (int d = 0; d < DIM; d++)
547                     {
548                         if (md->invMassPerDim[i][d] > 1.1*ALMOST_ZERO)
549                         {
550                             md->invMassPerDim[i][d] = md->invmass[i];
551                         }
552                     }
553                 }
554             }
555         }
556
557         /* Update the system mass for the change in lambda */
558         md->tmass  = L1*md->tmassA + lambda*md->tmassB;
559     }
560
561     md->lambda = lambda;
562 }