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