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