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