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