176ef3d751af858fb51c57c5e9a126e714bf2ef4
[alexxy/gromacs.git] / src / gromacs / mdlib / md_support.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) 2013,2014,2015,2016,2017 The GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #include "gmxpre.h"
40
41 #include "md_support.h"
42
43 #include <climits>
44 #include <cmath>
45
46 #include <algorithm>
47
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/coupling.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdlib/simulationsignal.h"
56 #include "gromacs/mdlib/stat.h"
57 #include "gromacs/mdlib/tgroup.h"
58 #include "gromacs/mdlib/update.h"
59 #include "gromacs/mdlib/vcm.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/df_history.h"
62 #include "gromacs/mdtypes/enerdata.h"
63 #include "gromacs/mdtypes/energyhistory.h"
64 #include "gromacs/mdtypes/forcerec.h"
65 #include "gromacs/mdtypes/group.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/pulling/pull.h"
72 #include "gromacs/timing/wallcycle.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/trajectory/trajectoryframe.h"
75 #include "gromacs/utility/arrayref.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.h"
79 #include "gromacs/utility/logger.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/snprintf.h"
82
83 static void calc_ke_part_normal(gmx::ArrayRef<const gmx::RVec> v,
84                                 const t_grpopts*               opts,
85                                 const t_mdatoms*               md,
86                                 gmx_ekindata_t*                ekind,
87                                 t_nrnb*                        nrnb,
88                                 gmx_bool                       bEkinAveVel)
89 {
90     int                         g;
91     gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
92
93     /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin.  Leap with AveVel is also
94        an option, but not supported now.
95        bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
96      */
97
98     // Now accumulate the partial global and groups ekin.
99     for (g = 0; (g < opts->ngtc); g++)
100     {
101         copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
102         if (bEkinAveVel)
103         {
104             clear_mat(tcstat[g].ekinf);
105             tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
106         }
107         else
108         {
109             clear_mat(tcstat[g].ekinh);
110         }
111     }
112     ekind->dekindl_old = ekind->dekindl;
113     int nthread        = gmx_omp_nthreads_get(emntUpdate);
114
115 #pragma omp parallel for num_threads(nthread) schedule(static)
116     for (int thread = 0; thread < nthread; thread++)
117     {
118         // This OpenMP only loops over arrays and does not call any functions
119         // or memory allocation. It should not be able to throw, so for now
120         // we do not need a try/catch wrapper.
121         int     start_t, end_t, n;
122         int     gt;
123         real    hm;
124         int     d, m;
125         matrix* ekin_sum;
126         real*   dekindl_sum;
127
128         start_t = ((thread + 0) * md->homenr) / nthread;
129         end_t   = ((thread + 1) * md->homenr) / nthread;
130
131         ekin_sum    = ekind->ekin_work[thread];
132         dekindl_sum = ekind->dekindl_work[thread];
133
134         for (gt = 0; gt < opts->ngtc; gt++)
135         {
136             clear_mat(ekin_sum[gt]);
137         }
138         *dekindl_sum = 0.0;
139
140         gt = 0;
141         for (n = start_t; n < end_t; n++)
142         {
143             if (md->cTC)
144             {
145                 gt = md->cTC[n];
146             }
147             hm = 0.5 * md->massT[n];
148
149             for (d = 0; (d < DIM); d++)
150             {
151                 for (m = 0; (m < DIM); m++)
152                 {
153                     /* if we're computing a full step velocity, v[d] has v(t).  Otherwise, v(t+dt/2) */
154                     ekin_sum[gt][m][d] += hm * v[n][m] * v[n][d];
155                 }
156             }
157             if (md->nMassPerturbed && md->bPerturbed[n])
158             {
159                 *dekindl_sum += 0.5 * (md->massB[n] - md->massA[n]) * iprod(v[n], v[n]);
160             }
161         }
162     }
163
164     ekind->dekindl = 0;
165     for (int thread = 0; thread < nthread; thread++)
166     {
167         for (g = 0; g < opts->ngtc; g++)
168         {
169             if (bEkinAveVel)
170             {
171                 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g], tcstat[g].ekinf);
172             }
173             else
174             {
175                 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g], tcstat[g].ekinh);
176             }
177         }
178
179         ekind->dekindl += *ekind->dekindl_work[thread];
180     }
181
182     inc_nrnb(nrnb, eNR_EKIN, md->homenr);
183 }
184
185 static void calc_ke_part_visc(const matrix                   box,
186                               gmx::ArrayRef<const gmx::RVec> x,
187                               gmx::ArrayRef<const gmx::RVec> v,
188                               const t_grpopts*               opts,
189                               const t_mdatoms*               md,
190                               gmx_ekindata_t*                ekind,
191                               t_nrnb*                        nrnb,
192                               gmx_bool                       bEkinAveVel)
193 {
194     int                         start = 0, homenr = md->homenr;
195     int                         g, d, n, m, gt = 0;
196     rvec                        v_corrt;
197     real                        hm;
198     gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
199     t_cos_acc*                  cosacc = &(ekind->cosacc);
200     real                        dekindl;
201     real                        fac, cosz;
202     double                      mvcos;
203
204     for (g = 0; g < opts->ngtc; g++)
205     {
206         copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
207         clear_mat(ekind->tcstat[g].ekinh);
208     }
209     ekind->dekindl_old = ekind->dekindl;
210
211     fac     = 2 * M_PI / box[ZZ][ZZ];
212     mvcos   = 0;
213     dekindl = 0;
214     for (n = start; n < start + homenr; n++)
215     {
216         if (md->cTC)
217         {
218             gt = md->cTC[n];
219         }
220         hm = 0.5 * md->massT[n];
221
222         /* Note that the times of x and v differ by half a step */
223         /* MRS -- would have to be changed for VV */
224         cosz = std::cos(fac * x[n][ZZ]);
225         /* Calculate the amplitude of the new velocity profile */
226         mvcos += 2 * cosz * md->massT[n] * v[n][XX];
227
228         copy_rvec(v[n], v_corrt);
229         /* Subtract the profile for the kinetic energy */
230         v_corrt[XX] -= cosz * cosacc->vcos;
231         for (d = 0; (d < DIM); d++)
232         {
233             for (m = 0; (m < DIM); m++)
234             {
235                 /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
236                 if (bEkinAveVel)
237                 {
238                     tcstat[gt].ekinf[m][d] += hm * v_corrt[m] * v_corrt[d];
239                 }
240                 else
241                 {
242                     tcstat[gt].ekinh[m][d] += hm * v_corrt[m] * v_corrt[d];
243                 }
244             }
245         }
246         if (md->nPerturbed && md->bPerturbed[n])
247         {
248             /* The minus sign here might be confusing.
249              * The kinetic contribution from dH/dl doesn't come from
250              * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
251              * where p are the momenta. The difference is only a minus sign.
252              */
253             dekindl -= 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
254         }
255     }
256     ekind->dekindl = dekindl;
257     cosacc->mvcos  = mvcos;
258
259     inc_nrnb(nrnb, eNR_EKIN, homenr);
260 }
261
262 static void calc_ke_part(gmx::ArrayRef<const gmx::RVec> x,
263                          gmx::ArrayRef<const gmx::RVec> v,
264                          const matrix                   box,
265                          const t_grpopts*               opts,
266                          const t_mdatoms*               md,
267                          gmx_ekindata_t*                ekind,
268                          t_nrnb*                        nrnb,
269                          gmx_bool                       bEkinAveVel)
270 {
271     if (ekind->cosacc.cos_accel == 0)
272     {
273         calc_ke_part_normal(v, opts, md, ekind, nrnb, bEkinAveVel);
274     }
275     else
276     {
277         calc_ke_part_visc(box, x, v, opts, md, ekind, nrnb, bEkinAveVel);
278     }
279 }
280
281 /* TODO Specialize this routine into init-time and loop-time versions?
282    e.g. bReadEkin is only true when restoring from checkpoint */
283 void compute_globals(gmx_global_stat*               gstat,
284                      t_commrec*                     cr,
285                      const t_inputrec*              ir,
286                      t_forcerec*                    fr,
287                      gmx_ekindata_t*                ekind,
288                      gmx::ArrayRef<const gmx::RVec> x,
289                      gmx::ArrayRef<const gmx::RVec> v,
290                      const matrix                   box,
291                      const t_mdatoms*               mdatoms,
292                      t_nrnb*                        nrnb,
293                      t_vcm*                         vcm,
294                      gmx_wallcycle_t                wcycle,
295                      gmx_enerdata_t*                enerd,
296                      tensor                         force_vir,
297                      tensor                         shake_vir,
298                      tensor                         total_vir,
299                      tensor                         pres,
300                      gmx::ArrayRef<real>            constraintsRmsdData,
301                      gmx::SimulationSignaller*      signalCoordinator,
302                      const matrix                   lastbox,
303                      int*                           totalNumberOfBondedInteractions,
304                      gmx_bool*                      bSumEkinhOld,
305                      const int                      flags)
306 {
307     gmx_bool bEner, bPres, bTemp;
308     gmx_bool bStopCM, bGStat, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
309     gmx_bool bCheckNumberOfBondedInteractions;
310     real     dvdl_ekin;
311
312     /* translate CGLO flags to gmx_booleans */
313     bStopCM                          = ((flags & CGLO_STOPCM) != 0);
314     bGStat                           = ((flags & CGLO_GSTAT) != 0);
315     bReadEkin                        = ((flags & CGLO_READEKIN) != 0);
316     bScaleEkin                       = ((flags & CGLO_SCALEEKIN) != 0);
317     bEner                            = ((flags & CGLO_ENERGY) != 0);
318     bTemp                            = ((flags & CGLO_TEMPERATURE) != 0);
319     bPres                            = ((flags & CGLO_PRESSURE) != 0);
320     bConstrain                       = ((flags & CGLO_CONSTRAINT) != 0);
321     bCheckNumberOfBondedInteractions = ((flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0);
322
323     /* we calculate a full state kinetic energy either with full-step velocity verlet
324        or half step where we need the pressure */
325
326     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
327
328     /* in initalization, it sums the shake virial in vv, and to
329        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
330
331     /* ########## Kinetic energy  ############## */
332
333     if (bTemp)
334     {
335         if (!bReadEkin)
336         {
337             calc_ke_part(x, v, box, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
338         }
339     }
340
341     /* Calculate center of mass velocity if necessary, also parallellized */
342     if (bStopCM)
343     {
344         calc_vcm_grp(*mdatoms, x, v, vcm);
345     }
346
347     if (bTemp || bStopCM || bPres || bEner || bConstrain || bCheckNumberOfBondedInteractions)
348     {
349         if (!bGStat)
350         {
351             /* We will not sum ekinh_old,
352              * so signal that we still have to do it.
353              */
354             *bSumEkinhOld = TRUE;
355         }
356         else
357         {
358             gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
359             if (PAR(cr))
360             {
361                 wallcycle_start(wcycle, ewcMoveE);
362                 global_stat(gstat,
363                             cr,
364                             enerd,
365                             force_vir,
366                             shake_vir,
367                             ir,
368                             ekind,
369                             constraintsRmsdData,
370                             bStopCM ? vcm : nullptr,
371                             signalBuffer.size(),
372                             signalBuffer.data(),
373                             totalNumberOfBondedInteractions,
374                             *bSumEkinhOld,
375                             flags);
376                 wallcycle_stop(wcycle, ewcMoveE);
377             }
378             signalCoordinator->finalizeSignals();
379             *bSumEkinhOld = FALSE;
380         }
381     }
382
383     if (bEner)
384     {
385         /* Calculate the amplitude of the cosine velocity profile */
386         ekind->cosacc.vcos = ekind->cosacc.mvcos / mdatoms->tmass;
387     }
388
389     if (bTemp)
390     {
391         /* Sum the kinetic energies of the groups & calc temp */
392         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
393         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
394            Leap with AveVel is not supported; it's not clear that it will actually work.
395            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
396            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
397          */
398         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin, bEkinAveVel, bScaleEkin);
399         enerd->dvdl_lin[efptMASS] = static_cast<double>(dvdl_ekin);
400
401         enerd->term[F_EKIN] = trace(ekind->ekin);
402     }
403
404     /* ########## Now pressure ############## */
405     // TODO: For the VV integrator bConstrain is needed in the conditional. This is confusing, so get rid of this.
406     if (bPres || bConstrain)
407     {
408         m_add(force_vir, shake_vir, total_vir);
409
410         /* Calculate pressure and apply LR correction if PPPM is used.
411          * Use the box from last timestep since we already called update().
412          */
413
414         enerd->term[F_PRES] = calc_pres(fr->pbcType, ir->nwall, lastbox, ekind->ekin, total_vir, pres);
415     }
416 }
417
418 static void min_zero(int* n, int i)
419 {
420     if (i > 0 && (*n == 0 || i < *n))
421     {
422         *n = i;
423     }
424 }
425
426 static int lcd4(int i1, int i2, int i3, int i4)
427 {
428     int nst;
429
430     nst = 0;
431     min_zero(&nst, i1);
432     min_zero(&nst, i2);
433     min_zero(&nst, i3);
434     min_zero(&nst, i4);
435     if (nst == 0)
436     {
437         gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
438     }
439
440     while (nst > 1
441            && ((i1 > 0 && i1 % nst != 0) || (i2 > 0 && i2 % nst != 0) || (i3 > 0 && i3 % nst != 0)
442                || (i4 > 0 && i4 % nst != 0)))
443     {
444         nst--;
445     }
446
447     return nst;
448 }
449
450 int computeGlobalCommunicationPeriod(const t_inputrec* ir)
451 {
452     int nstglobalcomm = 10;
453     {
454         // Set up the default behaviour
455         if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != TemperatureCoupling::No
456               || ir->epc != PressureCoupling::No))
457         {
458             /* The user didn't choose the period for anything
459                important, so we just make sure we can send signals and
460                write output suitably. */
461             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
462             {
463                 nstglobalcomm = ir->nstenergy;
464             }
465         }
466         else
467         {
468             /* The user has made a choice (perhaps implicitly), so we
469              * ensure that we do timely intra-simulation communication
470              * for (possibly) each of the four parts that care.
471              *
472              * TODO Does the Verlet scheme (+ DD) need any
473              * communication at nstlist steps? Is the use of nstlist
474              * here a leftover of the twin-range scheme? Can we remove
475              * nstlist when we remove the group scheme?
476              */
477             nstglobalcomm = lcd4(ir->nstcalcenergy,
478                                  ir->nstlist,
479                                  ir->etc != TemperatureCoupling::No ? ir->nsttcouple : 0,
480                                  ir->epc != PressureCoupling::No ? ir->nstpcouple : 0);
481         }
482     }
483     return nstglobalcomm;
484 }
485
486 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, const t_inputrec* ir, const t_commrec* cr)
487 {
488     const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
489
490     if (cr->nnodes > 1)
491     {
492         GMX_LOG(mdlog.info)
493                 .appendTextFormatted("Intra-simulation communication will occur every %d steps.\n",
494                                      nstglobalcomm);
495     }
496     return nstglobalcomm;
497 }
498
499 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
500 {
501     rvec *xp, *vp;
502
503     if (MASTER(cr) && *bLastStep)
504     {
505         fr->natoms = -1;
506     }
507     xp = fr->x;
508     vp = fr->v;
509     gmx_bcast(sizeof(*fr), fr, cr->mpi_comm_mygroup);
510     fr->x = xp;
511     fr->v = vp;
512
513     *bLastStep = (fr->natoms < 0);
514 }
515
516 // TODO Most of this logic seems to belong in the respective modules
517 void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
518 {
519     /* The entries in the state in the tpx file might not correspond
520      * with what is needed, so we correct this here.
521      */
522     state->flags = 0;
523     if (ir->efep != efepNO || ir->bExpanded)
524     {
525         state->flags |= (1 << estLAMBDA);
526         state->flags |= (1 << estFEPSTATE);
527     }
528     state->flags |= (1 << estX);
529     GMX_RELEASE_ASSERT(state->x.size() == state->natoms,
530                        "We should start a run with an initialized state->x");
531     if (EI_DYNAMICS(ir->eI))
532     {
533         state->flags |= (1 << estV);
534     }
535
536     state->nnhpres = 0;
537     if (ir->pbcType != PbcType::No)
538     {
539         state->flags |= (1 << estBOX);
540         if (inputrecPreserveShape(ir))
541         {
542             state->flags |= (1 << estBOX_REL);
543         }
544         if ((ir->epc == PressureCoupling::ParrinelloRahman) || (ir->epc == PressureCoupling::Mttk))
545         {
546             state->flags |= (1 << estBOXV);
547             if (!useModularSimulator)
548             {
549                 state->flags |= (1 << estPRES_PREV);
550             }
551         }
552         if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
553         {
554             state->nnhpres = 1;
555             state->flags |= (1 << estNHPRES_XI);
556             state->flags |= (1 << estNHPRES_VXI);
557             state->flags |= (1 << estSVIR_PREV);
558             state->flags |= (1 << estFVIR_PREV);
559             state->flags |= (1 << estVETA);
560             state->flags |= (1 << estVOL0);
561         }
562         if (ir->epc == PressureCoupling::Berendsen || ir->epc == PressureCoupling::CRescale)
563         {
564             state->flags |= (1 << estBAROS_INT);
565         }
566     }
567
568     if (ir->etc == TemperatureCoupling::NoseHoover)
569     {
570         state->flags |= (1 << estNH_XI);
571         state->flags |= (1 << estNH_VXI);
572     }
573
574     if (ir->etc == TemperatureCoupling::VRescale || ir->etc == TemperatureCoupling::Berendsen)
575     {
576         state->flags |= (1 << estTHERM_INT);
577     }
578
579     init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
580     init_ekinstate(&state->ekinstate, ir);
581
582     if (ir->bExpanded)
583     {
584         snew(state->dfhist, 1);
585         init_df_history(state->dfhist, ir->fepvals->n_lambda);
586     }
587
588     if (ir->pull && ir->pull->bSetPbcRefToPrevStepCOM)
589     {
590         state->flags |= (1 << estPULLCOMPREVSTEP);
591     }
592 }