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