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