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