Make SD stuff private for update.cpp
[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, 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/vec.h"
52 #include "gromacs/mdlib/coupling.h"
53 #include "gromacs/mdlib/dispersioncorrection.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/mdrunutility/multisim.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 // TODO move this to multi-sim module
85 bool multisim_int_all_are_equal(const gmx_multisim_t* ms, int64_t value)
86 {
87     bool     allValuesAreEqual = true;
88     int64_t* buf;
89
90     GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
91
92     snew(buf, ms->nsim);
93     /* send our value to all other master ranks, receive all of theirs */
94     buf[ms->sim] = value;
95     gmx_sumli_sim(ms->nsim, buf, ms);
96
97     for (int s = 0; s < ms->nsim; s++)
98     {
99         if (buf[s] != value)
100         {
101             allValuesAreEqual = false;
102             break;
103         }
104     }
105
106     sfree(buf);
107
108     return allValuesAreEqual;
109 }
110
111 int multisim_min(const gmx_multisim_t* ms, int nmin, int n)
112 {
113     int*     buf;
114     gmx_bool bPos, bEqual;
115     int      s, d;
116
117     snew(buf, ms->nsim);
118     buf[ms->sim] = n;
119     gmx_sumi_sim(ms->nsim, buf, ms);
120     bPos   = TRUE;
121     bEqual = TRUE;
122     for (s = 0; s < ms->nsim; s++)
123     {
124         bPos   = bPos && (buf[s] > 0);
125         bEqual = bEqual && (buf[s] == buf[0]);
126     }
127     if (bPos)
128     {
129         if (bEqual)
130         {
131             nmin = std::min(nmin, buf[0]);
132         }
133         else
134         {
135             /* Find the least common multiple */
136             for (d = 2; d < nmin; d++)
137             {
138                 s = 0;
139                 while (s < ms->nsim && d % buf[s] == 0)
140                 {
141                     s++;
142                 }
143                 if (s == ms->nsim)
144                 {
145                     /* We found the LCM and it is less than nmin */
146                     nmin = d;
147                     break;
148                 }
149             }
150         }
151     }
152     sfree(buf);
153
154     return nmin;
155 }
156
157 static void calc_ke_part_normal(gmx::ArrayRef<const gmx::RVec> v,
158                                 const t_grpopts*               opts,
159                                 const t_mdatoms*               md,
160                                 gmx_ekindata_t*                ekind,
161                                 t_nrnb*                        nrnb,
162                                 gmx_bool                       bEkinAveVel)
163 {
164     int                         g;
165     gmx::ArrayRef<t_grp_tcstat> tcstat  = ekind->tcstat;
166     gmx::ArrayRef<t_grp_acc>    grpstat = ekind->grpstat;
167
168     /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin.  Leap with AveVel is also
169        an option, but not supported now.
170        bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
171      */
172
173     /* group velocities are calculated in update_ekindata and
174      * accumulated in acumulate_groups.
175      * Now the partial global and groups ekin.
176      */
177     for (g = 0; (g < opts->ngtc); g++)
178     {
179         copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
180         if (bEkinAveVel)
181         {
182             clear_mat(tcstat[g].ekinf);
183             tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
184         }
185         else
186         {
187             clear_mat(tcstat[g].ekinh);
188         }
189     }
190     ekind->dekindl_old = ekind->dekindl;
191     int nthread        = gmx_omp_nthreads_get(emntUpdate);
192
193 #pragma omp parallel for num_threads(nthread) schedule(static)
194     for (int thread = 0; thread < nthread; thread++)
195     {
196         // This OpenMP only loops over arrays and does not call any functions
197         // or memory allocation. It should not be able to throw, so for now
198         // we do not need a try/catch wrapper.
199         int     start_t, end_t, n;
200         int     ga, gt;
201         rvec    v_corrt;
202         real    hm;
203         int     d, m;
204         matrix* ekin_sum;
205         real*   dekindl_sum;
206
207         start_t = ((thread + 0) * md->homenr) / nthread;
208         end_t   = ((thread + 1) * md->homenr) / nthread;
209
210         ekin_sum    = ekind->ekin_work[thread];
211         dekindl_sum = ekind->dekindl_work[thread];
212
213         for (gt = 0; gt < opts->ngtc; gt++)
214         {
215             clear_mat(ekin_sum[gt]);
216         }
217         *dekindl_sum = 0.0;
218
219         ga = 0;
220         gt = 0;
221         for (n = start_t; n < end_t; n++)
222         {
223             if (md->cACC)
224             {
225                 ga = md->cACC[n];
226             }
227             if (md->cTC)
228             {
229                 gt = md->cTC[n];
230             }
231             hm = 0.5 * md->massT[n];
232
233             for (d = 0; (d < DIM); d++)
234             {
235                 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
236             }
237             for (d = 0; (d < DIM); d++)
238             {
239                 for (m = 0; (m < DIM); m++)
240                 {
241                     /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
242                     ekin_sum[gt][m][d] += hm * v_corrt[m] * v_corrt[d];
243                 }
244             }
245             if (md->nMassPerturbed && md->bPerturbed[n])
246             {
247                 *dekindl_sum += 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
248             }
249         }
250     }
251
252     ekind->dekindl = 0;
253     for (int thread = 0; thread < nthread; thread++)
254     {
255         for (g = 0; g < opts->ngtc; g++)
256         {
257             if (bEkinAveVel)
258             {
259                 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g], tcstat[g].ekinf);
260             }
261             else
262             {
263                 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g], tcstat[g].ekinh);
264             }
265         }
266
267         ekind->dekindl += *ekind->dekindl_work[thread];
268     }
269
270     inc_nrnb(nrnb, eNR_EKIN, md->homenr);
271 }
272
273 static void calc_ke_part_visc(const matrix                   box,
274                               gmx::ArrayRef<const gmx::RVec> x,
275                               gmx::ArrayRef<const gmx::RVec> v,
276                               const t_grpopts*               opts,
277                               const t_mdatoms*               md,
278                               gmx_ekindata_t*                ekind,
279                               t_nrnb*                        nrnb,
280                               gmx_bool                       bEkinAveVel)
281 {
282     int                         start = 0, homenr = md->homenr;
283     int                         g, d, n, m, gt = 0;
284     rvec                        v_corrt;
285     real                        hm;
286     gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
287     t_cos_acc*                  cosacc = &(ekind->cosacc);
288     real                        dekindl;
289     real                        fac, cosz;
290     double                      mvcos;
291
292     for (g = 0; g < opts->ngtc; g++)
293     {
294         copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
295         clear_mat(ekind->tcstat[g].ekinh);
296     }
297     ekind->dekindl_old = ekind->dekindl;
298
299     fac     = 2 * M_PI / box[ZZ][ZZ];
300     mvcos   = 0;
301     dekindl = 0;
302     for (n = start; n < start + homenr; n++)
303     {
304         if (md->cTC)
305         {
306             gt = md->cTC[n];
307         }
308         hm = 0.5 * md->massT[n];
309
310         /* Note that the times of x and v differ by half a step */
311         /* MRS -- would have to be changed for VV */
312         cosz = std::cos(fac * x[n][ZZ]);
313         /* Calculate the amplitude of the new velocity profile */
314         mvcos += 2 * cosz * md->massT[n] * v[n][XX];
315
316         copy_rvec(v[n], v_corrt);
317         /* Subtract the profile for the kinetic energy */
318         v_corrt[XX] -= cosz * cosacc->vcos;
319         for (d = 0; (d < DIM); d++)
320         {
321             for (m = 0; (m < DIM); m++)
322             {
323                 /* if we're computing a full step velocity, v_corrt[d] has v(t).  Otherwise, v(t+dt/2) */
324                 if (bEkinAveVel)
325                 {
326                     tcstat[gt].ekinf[m][d] += hm * v_corrt[m] * v_corrt[d];
327                 }
328                 else
329                 {
330                     tcstat[gt].ekinh[m][d] += hm * v_corrt[m] * v_corrt[d];
331                 }
332             }
333         }
334         if (md->nPerturbed && md->bPerturbed[n])
335         {
336             /* The minus sign here might be confusing.
337              * The kinetic contribution from dH/dl doesn't come from
338              * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
339              * where p are the momenta. The difference is only a minus sign.
340              */
341             dekindl -= 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
342         }
343     }
344     ekind->dekindl = dekindl;
345     cosacc->mvcos  = mvcos;
346
347     inc_nrnb(nrnb, eNR_EKIN, homenr);
348 }
349
350 static void calc_ke_part(gmx::ArrayRef<const gmx::RVec> x,
351                          gmx::ArrayRef<const gmx::RVec> v,
352                          const matrix                   box,
353                          const t_grpopts*               opts,
354                          const t_mdatoms*               md,
355                          gmx_ekindata_t*                ekind,
356                          t_nrnb*                        nrnb,
357                          gmx_bool                       bEkinAveVel)
358 {
359     if (ekind->cosacc.cos_accel == 0)
360     {
361         calc_ke_part_normal(v, opts, md, ekind, nrnb, bEkinAveVel);
362     }
363     else
364     {
365         calc_ke_part_visc(box, x, v, opts, md, ekind, nrnb, bEkinAveVel);
366     }
367 }
368
369 /* TODO Specialize this routine into init-time and loop-time versions?
370    e.g. bReadEkin is only true when restoring from checkpoint */
371 void compute_globals(gmx_global_stat*               gstat,
372                      t_commrec*                     cr,
373                      const t_inputrec*              ir,
374                      t_forcerec*                    fr,
375                      gmx_ekindata_t*                ekind,
376                      gmx::ArrayRef<const gmx::RVec> x,
377                      gmx::ArrayRef<const gmx::RVec> v,
378                      const matrix                   box,
379                      real                           vdwLambda,
380                      const t_mdatoms*               mdatoms,
381                      t_nrnb*                        nrnb,
382                      t_vcm*                         vcm,
383                      gmx_wallcycle_t                wcycle,
384                      gmx_enerdata_t*                enerd,
385                      tensor                         force_vir,
386                      tensor                         shake_vir,
387                      tensor                         total_vir,
388                      tensor                         pres,
389                      gmx::Constraints*              constr,
390                      gmx::SimulationSignaller*      signalCoordinator,
391                      const matrix                   lastbox,
392                      int*                           totalNumberOfBondedInteractions,
393                      gmx_bool*                      bSumEkinhOld,
394                      const int                      flags)
395 {
396     gmx_bool bEner, bPres, bTemp;
397     gmx_bool bStopCM, bGStat, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
398     gmx_bool bCheckNumberOfBondedInteractions;
399     real     dvdl_ekin;
400
401     /* translate CGLO flags to gmx_booleans */
402     bStopCM                          = ((flags & CGLO_STOPCM) != 0);
403     bGStat                           = ((flags & CGLO_GSTAT) != 0);
404     bReadEkin                        = ((flags & CGLO_READEKIN) != 0);
405     bScaleEkin                       = ((flags & CGLO_SCALEEKIN) != 0);
406     bEner                            = ((flags & CGLO_ENERGY) != 0);
407     bTemp                            = ((flags & CGLO_TEMPERATURE) != 0);
408     bPres                            = ((flags & CGLO_PRESSURE) != 0);
409     bConstrain                       = ((flags & CGLO_CONSTRAINT) != 0);
410     bCheckNumberOfBondedInteractions = ((flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0);
411
412     /* we calculate a full state kinetic energy either with full-step velocity verlet
413        or half step where we need the pressure */
414
415     bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
416
417     /* in initalization, it sums the shake virial in vv, and to
418        sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
419
420     /* ########## Kinetic energy  ############## */
421
422     if (bTemp)
423     {
424         /* Non-equilibrium MD: this is parallellized, but only does communication
425          * when there really is NEMD.
426          */
427
428         if (PAR(cr) && (ekind->bNEMD))
429         {
430             accumulate_u(cr, &(ir->opts), ekind);
431         }
432         if (!bReadEkin)
433         {
434             calc_ke_part(x, v, box, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
435         }
436     }
437
438     /* Calculate center of mass velocity if necessary, also parallellized */
439     if (bStopCM)
440     {
441         calc_vcm_grp(*mdatoms, x, v, vcm);
442     }
443
444     if (bTemp || bStopCM || bPres || bEner || bConstrain || bCheckNumberOfBondedInteractions)
445     {
446         if (!bGStat)
447         {
448             /* We will not sum ekinh_old,
449              * so signal that we still have to do it.
450              */
451             *bSumEkinhOld = TRUE;
452         }
453         else
454         {
455             gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
456             if (PAR(cr))
457             {
458                 wallcycle_start(wcycle, ewcMoveE);
459                 global_stat(gstat, cr, enerd, force_vir, shake_vir, ir, ekind, constr,
460                             bStopCM ? vcm : nullptr, signalBuffer.size(), signalBuffer.data(),
461                             totalNumberOfBondedInteractions, *bSumEkinhOld, flags);
462                 wallcycle_stop(wcycle, ewcMoveE);
463             }
464             signalCoordinator->finalizeSignals();
465             *bSumEkinhOld = FALSE;
466         }
467     }
468
469     if (bEner)
470     {
471         /* Calculate the amplitude of the cosine velocity profile */
472         ekind->cosacc.vcos = ekind->cosacc.mvcos / mdatoms->tmass;
473     }
474
475     if (bTemp)
476     {
477         /* Sum the kinetic energies of the groups & calc temp */
478         /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
479         /* three maincase:  VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
480            Leap with AveVel is not supported; it's not clear that it will actually work.
481            bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
482            If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
483          */
484         enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin, bEkinAveVel, bScaleEkin);
485         enerd->dvdl_lin[efptMASS] = static_cast<double>(dvdl_ekin);
486
487         enerd->term[F_EKIN] = trace(ekind->ekin);
488
489         for (auto& dhdl : enerd->dhdlLambda)
490         {
491             dhdl += enerd->dvdl_lin[efptMASS];
492         }
493     }
494
495     /* ########## Now pressure ############## */
496     // TODO: For the VV integrator bConstrain is needed in the conditional. This is confusing, so get rid of this.
497     if (bPres || bConstrain)
498     {
499         m_add(force_vir, shake_vir, total_vir);
500
501         /* Calculate pressure and apply LR correction if PPPM is used.
502          * Use the box from last timestep since we already called update().
503          */
504
505         enerd->term[F_PRES] = calc_pres(fr->pbcType, ir->nwall, lastbox, ekind->ekin, total_vir, pres);
506     }
507
508     /* ##########  Long range energy information ###### */
509     if ((bEner || bPres) && fr->dispersionCorrection)
510     {
511         /* Calculate long range corrections to pressure and energy */
512         /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
513            and computes enerd->term[F_DISPCORR].  Also modifies the
514            total_vir and pres tensors */
515
516         const DispersionCorrection::Correction correction =
517                 fr->dispersionCorrection->calculate(lastbox, vdwLambda);
518
519         if (bEner)
520         {
521             enerd->term[F_DISPCORR] = correction.energy;
522             enerd->term[F_EPOT] += correction.energy;
523             enerd->term[F_DVDL_VDW] += correction.dvdl;
524         }
525         if (bPres)
526         {
527             correction.correctVirial(total_vir);
528             correction.correctPressure(pres);
529             enerd->term[F_PDISPCORR] = correction.pressure;
530             enerd->term[F_PRES] += correction.pressure;
531         }
532     }
533 }
534
535 void setCurrentLambdasRerun(int64_t           step,
536                             const t_lambda*   fepvals,
537                             const t_trxframe* rerun_fr,
538                             const double*     lam0,
539                             t_state*          globalState)
540 {
541     GMX_RELEASE_ASSERT(globalState != nullptr,
542                        "setCurrentLambdasGlobalRerun should be called with a valid state object");
543
544     if (rerun_fr->bLambda)
545     {
546         if (fepvals->delta_lambda == 0)
547         {
548             globalState->lambda[efptFEP] = rerun_fr->lambda;
549         }
550         else
551         {
552             /* find out between which two value of lambda we should be */
553             real frac      = step * fepvals->delta_lambda;
554             int  fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
555             /* interpolate between this state and the next */
556             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
557             frac = frac * fepvals->n_lambda - fep_state;
558             for (int i = 0; i < efptNR; i++)
559             {
560                 globalState->lambda[i] =
561                         lam0[i] + (fepvals->all_lambda[i][fep_state])
562                         + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
563             }
564         }
565     }
566     else if (rerun_fr->bFepState)
567     {
568         globalState->fep_state = rerun_fr->fep_state;
569         for (int i = 0; i < efptNR; i++)
570         {
571             globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
572         }
573     }
574 }
575
576 void setCurrentLambdasLocal(const int64_t       step,
577                             const t_lambda*     fepvals,
578                             const double*       lam0,
579                             gmx::ArrayRef<real> lambda,
580                             const int           currentFEPState)
581 /* find the current lambdas.  If rerunning, we either read in a state, or a lambda value,
582    requiring different logic. */
583 {
584     if (fepvals->delta_lambda != 0)
585     {
586         /* find out between which two value of lambda we should be */
587         real frac = step * fepvals->delta_lambda;
588         if (fepvals->n_lambda > 0)
589         {
590             int fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
591             /* interpolate between this state and the next */
592             /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
593             frac = frac * fepvals->n_lambda - fep_state;
594             for (int i = 0; i < efptNR; i++)
595             {
596                 lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state])
597                             + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
598             }
599         }
600         else
601         {
602             for (int i = 0; i < efptNR; i++)
603             {
604                 lambda[i] = lam0[i] + frac;
605             }
606         }
607     }
608     else
609     {
610         /* if < 0, fep_state was never defined, and we should not set lambda from the state */
611         if (currentFEPState > -1)
612         {
613             for (int i = 0; i < efptNR; i++)
614             {
615                 lambda[i] = fepvals->all_lambda[i][currentFEPState];
616             }
617         }
618     }
619 }
620
621 static void min_zero(int* n, int i)
622 {
623     if (i > 0 && (*n == 0 || i < *n))
624     {
625         *n = i;
626     }
627 }
628
629 static int lcd4(int i1, int i2, int i3, int i4)
630 {
631     int nst;
632
633     nst = 0;
634     min_zero(&nst, i1);
635     min_zero(&nst, i2);
636     min_zero(&nst, i3);
637     min_zero(&nst, i4);
638     if (nst == 0)
639     {
640         gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
641     }
642
643     while (nst > 1
644            && ((i1 > 0 && i1 % nst != 0) || (i2 > 0 && i2 % nst != 0) || (i3 > 0 && i3 % nst != 0)
645                || (i4 > 0 && i4 % nst != 0)))
646     {
647         nst--;
648     }
649
650     return nst;
651 }
652
653 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr)
654 {
655     int nstglobalcomm;
656     {
657         // Set up the default behaviour
658         if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != etcNO || ir->epc != epcNO))
659         {
660             /* The user didn't choose the period for anything
661                important, so we just make sure we can send signals and
662                write output suitably. */
663             nstglobalcomm = 10;
664             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
665             {
666                 nstglobalcomm = ir->nstenergy;
667             }
668         }
669         else
670         {
671             /* The user has made a choice (perhaps implicitly), so we
672              * ensure that we do timely intra-simulation communication
673              * for (possibly) each of the four parts that care.
674              *
675              * TODO Does the Verlet scheme (+ DD) need any
676              * communication at nstlist steps? Is the use of nstlist
677              * here a leftover of the twin-range scheme? Can we remove
678              * nstlist when we remove the group scheme?
679              */
680             nstglobalcomm = lcd4(ir->nstcalcenergy, ir->nstlist, ir->etc != etcNO ? ir->nsttcouple : 0,
681                                  ir->epc != epcNO ? ir->nstpcouple : 0);
682         }
683     }
684
685     // TODO change this behaviour. Instead grompp should print
686     // a (performance) note and mdrun should not change ir.
687     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
688     {
689         GMX_LOG(mdlog.warning)
690                 .asParagraph()
691                 .appendTextFormatted("WARNING: Changing nstcomm from %d to %d", ir->nstcomm, nstglobalcomm);
692         ir->nstcomm = nstglobalcomm;
693     }
694
695     if (cr->nnodes > 1)
696     {
697         GMX_LOG(mdlog.info)
698                 .appendTextFormatted("Intra-simulation communication will occur every %d steps.\n",
699                                      nstglobalcomm);
700     }
701     return nstglobalcomm;
702 }
703
704 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
705 {
706     rvec *xp, *vp;
707
708     if (MASTER(cr) && *bLastStep)
709     {
710         fr->natoms = -1;
711     }
712     xp = fr->x;
713     vp = fr->v;
714     gmx_bcast(sizeof(*fr), fr, cr->mpi_comm_mygroup);
715     fr->x = xp;
716     fr->v = vp;
717
718     *bLastStep = (fr->natoms < 0);
719 }
720
721 // TODO Most of this logic seems to belong in the respective modules
722 void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
723 {
724     /* The entries in the state in the tpx file might not correspond
725      * with what is needed, so we correct this here.
726      */
727     state->flags = 0;
728     if (ir->efep != efepNO || ir->bExpanded)
729     {
730         state->flags |= (1 << estLAMBDA);
731         state->flags |= (1 << estFEPSTATE);
732     }
733     state->flags |= (1 << estX);
734     GMX_RELEASE_ASSERT(state->x.size() == state->natoms,
735                        "We should start a run with an initialized state->x");
736     if (EI_DYNAMICS(ir->eI))
737     {
738         state->flags |= (1 << estV);
739     }
740
741     state->nnhpres = 0;
742     if (ir->pbcType != PbcType::No)
743     {
744         state->flags |= (1 << estBOX);
745         if (inputrecPreserveShape(ir))
746         {
747             state->flags |= (1 << estBOX_REL);
748         }
749         if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
750         {
751             state->flags |= (1 << estBOXV);
752             if (!useModularSimulator)
753             {
754                 state->flags |= (1 << estPRES_PREV);
755             }
756         }
757         if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
758         {
759             state->nnhpres = 1;
760             state->flags |= (1 << estNHPRES_XI);
761             state->flags |= (1 << estNHPRES_VXI);
762             state->flags |= (1 << estSVIR_PREV);
763             state->flags |= (1 << estFVIR_PREV);
764             state->flags |= (1 << estVETA);
765             state->flags |= (1 << estVOL0);
766         }
767         if (ir->epc == epcBERENDSEN)
768         {
769             state->flags |= (1 << estBAROS_INT);
770         }
771     }
772
773     if (ir->etc == etcNOSEHOOVER)
774     {
775         state->flags |= (1 << estNH_XI);
776         state->flags |= (1 << estNH_VXI);
777     }
778
779     if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
780     {
781         state->flags |= (1 << estTHERM_INT);
782     }
783
784     init_gtc_state(state, state->ngtc, state->nnhpres,
785                    ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
786     init_ekinstate(&state->ekinstate, ir);
787
788     if (ir->bExpanded)
789     {
790         snew(state->dfhist, 1);
791         init_df_history(state->dfhist, ir->fepvals->n_lambda);
792     }
793
794     if (ir->pull && ir->pull->bSetPbcRefToPrevStepCOM)
795     {
796         state->flags |= (1 << estPULLCOMPREVSTEP);
797     }
798 }