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