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