Use more const references in topology code
[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_.resizeWithPadding(newSize);
105     mdatoms_->chargeA = chargeA_.data();
106 }
107
108 void MDAtoms::reserve(int newCapacity)
109 {
110     chargeA_.reserveWithPadding(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 (getGroupType(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) != 0);
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     int                   nthreads gmx_unused;
211
212     bLJPME = EVDW_PME(ir->vdwtype);
213
214     opts = &ir->opts;
215
216     const gmx_groups_t &groups = mtop->groups;
217
218     auto                md = mdAtoms->mdatoms();
219     /* nindex>=0 indicates DD where we use an index */
220     if (nindex >= 0)
221     {
222         md->nr = nindex;
223     }
224     else
225     {
226         md->nr = mtop->natoms;
227     }
228
229     if (md->nr > md->nalloc)
230     {
231         md->nalloc = over_alloc_dd(md->nr);
232
233         if (md->nMassPerturbed)
234         {
235             srenew(md->massA, md->nalloc);
236             srenew(md->massB, md->nalloc);
237         }
238         srenew(md->massT, md->nalloc);
239         /* The SIMD version of the integrator needs this aligned and padded.
240          * The padding needs to be with zeros, which we set later below.
241          */
242         gmx::AlignedAllocationPolicy::free(md->invmass);
243         md->invmass = new(gmx::AlignedAllocationPolicy::malloc((md->nalloc + GMX_REAL_MAX_SIMD_WIDTH)*sizeof(*md->invmass)))real;
244         srenew(md->invMassPerDim, md->nalloc);
245         // TODO eventually we will have vectors and just resize
246         // everything, but for now the semantics of md->nalloc being
247         // the capacity are preserved by keeping vectors within
248         // mdAtoms having the same properties as the other arrays.
249         mdAtoms->reserve(md->nalloc);
250         mdAtoms->resize(md->nr);
251         srenew(md->typeA, md->nalloc);
252         if (md->nPerturbed)
253         {
254             srenew(md->chargeB, md->nalloc);
255             srenew(md->typeB, md->nalloc);
256         }
257         if (bLJPME)
258         {
259             srenew(md->sqrt_c6A, md->nalloc);
260             srenew(md->sigmaA, md->nalloc);
261             srenew(md->sigma3A, md->nalloc);
262             if (md->nPerturbed)
263             {
264                 srenew(md->sqrt_c6B, md->nalloc);
265                 srenew(md->sigmaB, md->nalloc);
266                 srenew(md->sigma3B, md->nalloc);
267             }
268         }
269         srenew(md->ptype, md->nalloc);
270         if (opts->ngtc > 1)
271         {
272             srenew(md->cTC, md->nalloc);
273             /* We always copy cTC with domain decomposition */
274         }
275         srenew(md->cENER, md->nalloc);
276         if (opts->ngacc > 1)
277         {
278             srenew(md->cACC, md->nalloc);
279         }
280         if (opts->nFreeze &&
281             (opts->ngfrz > 1 ||
282              opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
283         {
284             srenew(md->cFREEZE, md->nalloc);
285         }
286         if (md->bVCMgrps)
287         {
288             srenew(md->cVCM, md->nalloc);
289         }
290         if (md->bOrires)
291         {
292             srenew(md->cORF, md->nalloc);
293         }
294         if (md->nPerturbed)
295         {
296             srenew(md->bPerturbed, md->nalloc);
297         }
298
299         /* Note that these user t_mdatoms array pointers are NULL
300          * when there is only one group present.
301          * Therefore, when adding code, the user should use something like:
302          * gprnrU1 = (md->cU1==NULL ? 0 : md->cU1[localatindex])
303          */
304         if (mtop->groups.grpnr[egcUser1] != nullptr)
305         {
306             srenew(md->cU1, md->nalloc);
307         }
308         if (mtop->groups.grpnr[egcUser2] != nullptr)
309         {
310             srenew(md->cU2, md->nalloc);
311         }
312
313         if (ir->bQMMM)
314         {
315             srenew(md->bQM, md->nalloc);
316         }
317     }
318
319     int molb = 0;
320
321     nthreads = gmx_omp_nthreads_get(emntDefault);
322 #pragma omp parallel for num_threads(nthreads) schedule(static) firstprivate(molb)
323     for (int i = 0; i < md->nr; i++)
324     {
325         try
326         {
327             int      g, ag;
328             real     mA, mB, fac;
329             real     c6, c12;
330
331             if (index == nullptr)
332             {
333                 ag = i;
334             }
335             else
336             {
337                 ag = index[i];
338             }
339             const t_atom &atom = mtopGetAtomParameters(mtop, ag, &molb);
340
341             if (md->cFREEZE)
342             {
343                 md->cFREEZE[i] = getGroupType(groups, egcFREEZE, ag);
344             }
345             if (EI_ENERGY_MINIMIZATION(ir->eI))
346             {
347                 /* Displacement is proportional to F, masses used for constraints */
348                 mA = 1.0;
349                 mB = 1.0;
350             }
351             else if (ir->eI == eiBD)
352             {
353                 /* With BD the physical masses are irrelevant.
354                  * To keep the code simple we use most of the normal MD code path
355                  * for BD. Thus for constraining the masses should be proportional
356                  * to the friction coefficient. We set the absolute value such that
357                  * m/2<(dx/dt)^2> = m/2*2kT/fric*dt = kT/2 => m=fric*dt/2
358                  * Then if we set the (meaningless) velocity to v=dx/dt, we get the
359                  * correct kinetic energy and temperature using the usual code path.
360                  * Thus with BD v*dt will give the displacement and the reported
361                  * temperature can signal bad integration (too large time step).
362                  */
363                 if (ir->bd_fric > 0)
364                 {
365                     mA = 0.5*ir->bd_fric*ir->delta_t;
366                     mB = 0.5*ir->bd_fric*ir->delta_t;
367                 }
368                 else
369                 {
370                     /* The friction coefficient is mass/tau_t */
371                     fac = ir->delta_t/opts->tau_t[md->cTC ? groups.grpnr[egcTC][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.grpnr[egcTC][ag];
469             }
470             md->cENER[i]    = getGroupType(groups, egcENER, ag);
471             if (md->cACC)
472             {
473                 md->cACC[i]   = groups.grpnr[egcACC][ag];
474             }
475             if (md->cVCM)
476             {
477                 md->cVCM[i]       = groups.grpnr[egcVCM][ag];
478             }
479             if (md->cORF)
480             {
481                 md->cORF[i]       = getGroupType(groups, egcORFIT, ag);
482             }
483
484             if (md->cU1)
485             {
486                 md->cU1[i]        = groups.grpnr[egcUser1][ag];
487             }
488             if (md->cU2)
489             {
490                 md->cU2[i]        = groups.grpnr[egcUser2][ag];
491             }
492
493             if (ir->bQMMM)
494             {
495                 if (groups.grpnr[egcQMMM] == nullptr ||
496                     groups.grpnr[egcQMMM][ag] < groups.grps[egcQMMM].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 }