Create separate module for trajectory data
[gromacs.git] / src / gromacs / mdlib / md_support.cpp
blob539ee1cba2896d664ec6dd85c0258a2c12de272c
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "md_support.h"
42 #include <algorithm>
44 #include "gromacs/domdec/domdec.h"
45 #include "gromacs/gmxlib/md_logging.h"
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/gmxlib/nrnb.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/mdrun.h"
50 #include "gromacs/mdlib/mdrun_signalling.h"
51 #include "gromacs/mdlib/tgroup.h"
52 #include "gromacs/mdlib/vcm.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/mdtypes/df_history.h"
55 #include "gromacs/mdtypes/energyhistory.h"
56 #include "gromacs/mdtypes/group.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/state.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/timing/wallcycle.h"
62 #include "gromacs/topology/mtop_util.h"
63 #include "gromacs/trajectory/trajectoryframe.h"
64 #include "gromacs/utility/arrayref.h"
65 #include "gromacs/utility/cstringutil.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/smalloc.h"
68 #include "gromacs/utility/snprintf.h"
70 /* check which of the multisim simulations has the shortest number of
71 steps and return that number of nsteps */
72 gmx_int64_t get_multisim_nsteps(const t_commrec *cr,
73 gmx_int64_t nsteps)
75 gmx_int64_t steps_out;
77 if (MASTER(cr))
79 gmx_int64_t *buf;
80 int s;
82 snew(buf, cr->ms->nsim);
84 buf[cr->ms->sim] = nsteps;
85 gmx_sumli_sim(cr->ms->nsim, buf, cr->ms);
87 steps_out = -1;
88 for (s = 0; s < cr->ms->nsim; s++)
90 /* find the smallest positive number */
91 if (buf[s] >= 0 && ((steps_out < 0) || (buf[s] < steps_out)) )
93 steps_out = buf[s];
96 sfree(buf);
98 /* if we're the limiting simulation, don't do anything */
99 if (steps_out >= 0 && steps_out < nsteps)
101 char strbuf[255];
102 snprintf(strbuf, 255, "Will stop simulation %%d after %s steps (another simulation will end then).\n", "%" GMX_PRId64);
103 fprintf(stderr, strbuf, cr->ms->sim, steps_out);
106 /* broadcast to non-masters */
107 gmx_bcast(sizeof(gmx_int64_t), &steps_out, cr);
108 return steps_out;
111 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
113 int *buf;
114 gmx_bool bPos, bEqual;
115 int s, d;
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++)
124 bPos = bPos && (buf[s] > 0);
125 bEqual = bEqual && (buf[s] == buf[0]);
127 if (bPos)
129 if (bEqual)
131 nmin = std::min(nmin, buf[0]);
133 else
135 /* Find the least common multiple */
136 for (d = 2; d < nmin; d++)
138 s = 0;
139 while (s < ms->nsim && d % buf[s] == 0)
141 s++;
143 if (s == ms->nsim)
145 /* We found the LCM and it is less than nmin */
146 nmin = d;
147 break;
152 sfree(buf);
154 return nmin;
157 int multisim_nstsimsync(const t_commrec *cr,
158 const t_inputrec *ir, int repl_ex_nst)
160 int nmin;
162 if (MASTER(cr))
164 nmin = INT_MAX;
165 nmin = multisim_min(cr->ms, nmin, ir->nstlist);
166 nmin = multisim_min(cr->ms, nmin, ir->nstcalcenergy);
167 nmin = multisim_min(cr->ms, nmin, repl_ex_nst);
168 if (nmin == INT_MAX)
170 gmx_fatal(FARGS, "Can not find an appropriate interval for inter-simulation communication, since nstlist, nstcalcenergy and -replex are all <= 0");
172 /* Avoid inter-simulation communication at every (second) step */
173 if (nmin <= 2)
175 nmin = 10;
179 gmx_bcast(sizeof(int), &nmin, cr);
181 return nmin;
184 void copy_coupling_state(t_state *statea, t_state *stateb,
185 gmx_ekindata_t *ekinda, gmx_ekindata_t *ekindb, t_grpopts* opts)
188 /* MRS note -- might be able to get rid of some of the arguments. Look over it when it's all debugged */
190 int i, j, nc;
192 /* Make sure we have enough space for x and v */
193 if (statea->nalloc > stateb->nalloc)
195 stateb->nalloc = statea->nalloc;
196 srenew(stateb->x, stateb->nalloc);
197 srenew(stateb->v, stateb->nalloc);
200 stateb->natoms = statea->natoms;
201 stateb->ngtc = statea->ngtc;
202 stateb->nnhpres = statea->nnhpres;
203 stateb->veta = statea->veta;
204 if (ekinda)
206 copy_mat(ekinda->ekin, ekindb->ekin);
207 for (i = 0; i < stateb->ngtc; i++)
209 ekindb->tcstat[i].T = ekinda->tcstat[i].T;
210 ekindb->tcstat[i].Th = ekinda->tcstat[i].Th;
211 copy_mat(ekinda->tcstat[i].ekinh, ekindb->tcstat[i].ekinh);
212 copy_mat(ekinda->tcstat[i].ekinf, ekindb->tcstat[i].ekinf);
213 ekindb->tcstat[i].ekinscalef_nhc = ekinda->tcstat[i].ekinscalef_nhc;
214 ekindb->tcstat[i].ekinscaleh_nhc = ekinda->tcstat[i].ekinscaleh_nhc;
215 ekindb->tcstat[i].vscale_nhc = ekinda->tcstat[i].vscale_nhc;
218 copy_rvecn(statea->x, stateb->x, 0, stateb->natoms);
219 copy_rvecn(statea->v, stateb->v, 0, stateb->natoms);
220 copy_mat(statea->box, stateb->box);
221 copy_mat(statea->box_rel, stateb->box_rel);
222 copy_mat(statea->boxv, stateb->boxv);
224 for (i = 0; i < stateb->ngtc; i++)
226 nc = i*opts->nhchainlength;
227 for (j = 0; j < opts->nhchainlength; j++)
229 stateb->nosehoover_xi[nc+j] = statea->nosehoover_xi[nc+j];
230 stateb->nosehoover_vxi[nc+j] = statea->nosehoover_vxi[nc+j];
233 if (stateb->nhpres_xi != NULL)
235 for (i = 0; i < stateb->nnhpres; i++)
237 nc = i*opts->nhchainlength;
238 for (j = 0; j < opts->nhchainlength; j++)
240 stateb->nhpres_xi[nc+j] = statea->nhpres_xi[nc+j];
241 stateb->nhpres_vxi[nc+j] = statea->nhpres_vxi[nc+j];
247 real compute_conserved_from_auxiliary(t_inputrec *ir, t_state *state, t_extmass *MassQ)
249 real quantity = 0;
250 switch (ir->etc)
252 case etcNO:
253 break;
254 case etcBERENDSEN:
255 break;
256 case etcNOSEHOOVER:
257 quantity = NPT_energy(ir, state, MassQ);
258 break;
259 case etcVRESCALE:
260 quantity = vrescale_energy(&(ir->opts), state->therm_integral);
261 break;
262 default:
263 break;
265 return quantity;
268 /* TODO Specialize this routine into init-time and loop-time versions?
269 e.g. bReadEkin is only true when restoring from checkpoint */
270 void compute_globals(FILE *fplog, gmx_global_stat_t gstat, t_commrec *cr, t_inputrec *ir,
271 t_forcerec *fr, gmx_ekindata_t *ekind,
272 t_state *state, t_mdatoms *mdatoms,
273 t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
274 gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
275 tensor pres, rvec mu_tot, gmx_constr_t constr,
276 struct gmx_signalling_t *gs, gmx_bool bInterSimGS,
277 matrix box, gmx_mtop_t *top_global,
278 gmx_bool *bSumEkinhOld, int flags)
280 tensor corr_vir, corr_pres;
281 gmx_bool bEner, bPres, bTemp;
282 gmx_bool bStopCM, bGStat,
283 bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
284 real prescorr, enercorr, dvdlcorr, dvdl_ekin;
286 /* translate CGLO flags to gmx_booleans */
287 bStopCM = flags & CGLO_STOPCM;
288 bGStat = flags & CGLO_GSTAT;
289 bReadEkin = (flags & CGLO_READEKIN);
290 bScaleEkin = (flags & CGLO_SCALEEKIN);
291 bEner = flags & CGLO_ENERGY;
292 bTemp = flags & CGLO_TEMPERATURE;
293 bPres = (flags & CGLO_PRESSURE);
294 bConstrain = (flags & CGLO_CONSTRAINT);
296 /* we calculate a full state kinetic energy either with full-step velocity verlet
297 or half step where we need the pressure */
299 bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
301 /* in initalization, it sums the shake virial in vv, and to
302 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
304 /* ########## Kinetic energy ############## */
306 if (bTemp)
308 /* Non-equilibrium MD: this is parallellized, but only does communication
309 * when there really is NEMD.
312 if (PAR(cr) && (ekind->bNEMD))
314 accumulate_u(cr, &(ir->opts), ekind);
316 if (!bReadEkin)
318 calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
322 /* Calculate center of mass velocity if necessary, also parallellized */
323 if (bStopCM)
325 calc_vcm_grp(0, mdatoms->homenr, mdatoms,
326 state->x, state->v, vcm);
329 if (bTemp || bStopCM || bPres || bEner || bConstrain)
331 if (!bGStat)
333 /* We will not sum ekinh_old,
334 * so signal that we still have to do it.
336 *bSumEkinhOld = TRUE;
339 else
341 gmx::ArrayRef<real> signalBuffer = prepareSignalBuffer(gs);
342 if (PAR(cr))
344 wallcycle_start(wcycle, ewcMoveE);
345 global_stat(fplog, gstat, cr, enerd, force_vir, shake_vir, mu_tot,
346 ir, ekind, constr, bStopCM ? vcm : NULL,
347 signalBuffer.size(), signalBuffer.data(),
348 top_global, state,
349 *bSumEkinhOld, flags);
350 wallcycle_stop(wcycle, ewcMoveE);
352 handleSignals(gs, cr, bInterSimGS);
353 *bSumEkinhOld = FALSE;
357 if (!ekind->bNEMD && debug && bTemp && (vcm->nr > 0))
359 correct_ekin(debug,
360 0, mdatoms->homenr,
361 state->v, vcm->group_p[0],
362 mdatoms->massT, mdatoms->tmass, ekind->ekin);
365 /* Do center of mass motion removal */
366 if (bStopCM)
368 check_cm_grp(fplog, vcm, ir, 1);
369 do_stopcm_grp(0, mdatoms->homenr, mdatoms->cVCM,
370 state->x, state->v, vcm);
371 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
374 if (bEner)
376 /* Calculate the amplitude of the cosine velocity profile */
377 ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
380 if (bTemp)
382 /* Sum the kinetic energies of the groups & calc temp */
383 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
384 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
385 Leap with AveVel is not supported; it's not clear that it will actually work.
386 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
387 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
389 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
390 bEkinAveVel, bScaleEkin);
391 enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
393 enerd->term[F_EKIN] = trace(ekind->ekin);
396 /* ########## Long range energy information ###### */
398 if (bEner || bPres || bConstrain)
400 calc_dispcorr(ir, fr, top_global->natoms, box, state->lambda[efptVDW],
401 corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
404 if (bEner)
406 enerd->term[F_DISPCORR] = enercorr;
407 enerd->term[F_EPOT] += enercorr;
408 enerd->term[F_DVDL_VDW] += dvdlcorr;
411 /* ########## Now pressure ############## */
412 if (bPres || bConstrain)
415 m_add(force_vir, shake_vir, total_vir);
417 /* Calculate pressure and apply LR correction if PPPM is used.
418 * Use the box from last timestep since we already called update().
421 enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
423 /* Calculate long range corrections to pressure and energy */
424 /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
425 and computes enerd->term[F_DISPCORR]. Also modifies the
426 total_vir and pres tesors */
428 m_add(total_vir, corr_vir, total_vir);
429 m_add(pres, corr_pres, pres);
430 enerd->term[F_PDISPCORR] = prescorr;
431 enerd->term[F_PRES] += prescorr;
435 void check_nst_param(FILE *fplog, t_commrec *cr,
436 const char *desc_nst, int nst,
437 const char *desc_p, int *p)
439 if (*p > 0 && *p % nst != 0)
441 /* Round up to the next multiple of nst */
442 *p = ((*p)/nst + 1)*nst;
443 md_print_warn(cr, fplog,
444 "NOTE: %s changes %s to %d\n", desc_nst, desc_p, *p);
448 void set_current_lambdas(gmx_int64_t step, t_lambda *fepvals, gmx_bool bRerunMD,
449 t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[])
450 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
451 requiring different logic. */
453 real frac;
454 int i, fep_state = 0;
455 if (bRerunMD)
457 if (rerun_fr->bLambda)
459 if (fepvals->delta_lambda == 0)
461 state_global->lambda[efptFEP] = rerun_fr->lambda;
462 for (i = 0; i < efptNR; i++)
464 if (i != efptFEP)
466 state->lambda[i] = state_global->lambda[i];
470 else
472 /* find out between which two value of lambda we should be */
473 frac = (step*fepvals->delta_lambda);
474 fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
475 /* interpolate between this state and the next */
476 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
477 frac = (frac*fepvals->n_lambda)-fep_state;
478 for (i = 0; i < efptNR; i++)
480 state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
481 frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
485 else if (rerun_fr->bFepState)
487 state_global->fep_state = rerun_fr->fep_state;
488 for (i = 0; i < efptNR; i++)
490 state_global->lambda[i] = fepvals->all_lambda[i][fep_state];
494 else
496 if (fepvals->delta_lambda != 0)
498 /* find out between which two value of lambda we should be */
499 frac = (step*fepvals->delta_lambda);
500 if (fepvals->n_lambda > 0)
502 fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
503 /* interpolate between this state and the next */
504 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
505 frac = (frac*fepvals->n_lambda)-fep_state;
506 for (i = 0; i < efptNR; i++)
508 state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
509 frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
512 else
514 for (i = 0; i < efptNR; i++)
516 state_global->lambda[i] = lam0[i] + frac;
520 else
522 if (state->fep_state > 0)
524 state_global->fep_state = state->fep_state; /* state->fep is the one updated by bExpanded */
525 for (i = 0; i < efptNR; i++)
527 state_global->lambda[i] = fepvals->all_lambda[i][state_global->fep_state];
532 for (i = 0; i < efptNR; i++)
534 state->lambda[i] = state_global->lambda[i];
538 static void min_zero(int *n, int i)
540 if (i > 0 && (*n == 0 || i < *n))
542 *n = i;
546 static int lcd4(int i1, int i2, int i3, int i4)
548 int nst;
550 nst = 0;
551 min_zero(&nst, i1);
552 min_zero(&nst, i2);
553 min_zero(&nst, i3);
554 min_zero(&nst, i4);
555 if (nst == 0)
557 gmx_incons("All 4 inputs for determininig nstglobalcomm are <= 0");
560 while (nst > 1 && ((i1 > 0 && i1 % nst != 0) ||
561 (i2 > 0 && i2 % nst != 0) ||
562 (i3 > 0 && i3 % nst != 0) ||
563 (i4 > 0 && i4 % nst != 0)))
565 nst--;
568 return nst;
571 int check_nstglobalcomm(FILE *fplog, t_commrec *cr,
572 int nstglobalcomm, t_inputrec *ir)
574 if (!EI_DYNAMICS(ir->eI))
576 nstglobalcomm = 1;
579 if (nstglobalcomm == -1)
581 if (!(ir->nstcalcenergy > 0 ||
582 ir->nstlist > 0 ||
583 ir->etc != etcNO ||
584 ir->epc != epcNO))
586 nstglobalcomm = 10;
587 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
589 nstglobalcomm = ir->nstenergy;
592 else
594 /* Ensure that we do timely global communication for
595 * (possibly) each of the four following options.
597 nstglobalcomm = lcd4(ir->nstcalcenergy,
598 ir->nstlist,
599 ir->etc != etcNO ? ir->nsttcouple : 0,
600 ir->epc != epcNO ? ir->nstpcouple : 0);
603 else
605 if (ir->nstlist > 0 &&
606 nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
608 nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
609 md_print_warn(cr, fplog, "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d\n", nstglobalcomm);
611 if (ir->nstcalcenergy > 0)
613 check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
614 "nstcalcenergy", &ir->nstcalcenergy);
616 if (ir->etc != etcNO && ir->nsttcouple > 0)
618 check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
619 "nsttcouple", &ir->nsttcouple);
621 if (ir->epc != epcNO && ir->nstpcouple > 0)
623 check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
624 "nstpcouple", &ir->nstpcouple);
627 check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
628 "nstenergy", &ir->nstenergy);
630 check_nst_param(fplog, cr, "-gcom", nstglobalcomm,
631 "nstlog", &ir->nstlog);
634 if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
636 md_print_warn(cr, fplog, "WARNING: Changing nstcomm from %d to %d\n",
637 ir->nstcomm, nstglobalcomm);
638 ir->nstcomm = nstglobalcomm;
641 return nstglobalcomm;
644 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
645 gmx_bool *bNotLastFrame)
647 rvec *xp, *vp;
649 if (MASTER(cr) && !*bNotLastFrame)
651 fr->natoms = -1;
653 xp = fr->x;
654 vp = fr->v;
655 gmx_bcast(sizeof(*fr), fr, cr);
656 fr->x = xp;
657 fr->v = vp;
659 *bNotLastFrame = (fr->natoms >= 0);
663 void set_state_entries(t_state *state, const t_inputrec *ir)
665 /* The entries in the state in the tpx file might not correspond
666 * with what is needed, so we correct this here.
668 state->flags = 0;
669 if (ir->efep != efepNO || ir->bExpanded)
671 state->flags |= (1<<estLAMBDA);
672 state->flags |= (1<<estFEPSTATE);
674 state->flags |= (1<<estX);
675 if (state->lambda == NULL)
677 snew(state->lambda, efptNR);
679 if (state->x == NULL)
681 snew(state->x, state->nalloc);
683 if (EI_DYNAMICS(ir->eI))
685 state->flags |= (1<<estV);
686 if (state->v == NULL)
688 snew(state->v, state->nalloc);
691 if (ir->eI == eiSD2)
693 state->flags |= (1<<estSDX);
694 if (state->sd_X == NULL)
696 /* sd_X is not stored in the tpx file, so we need to allocate it */
697 snew(state->sd_X, state->nalloc);
700 if (ir->eI == eiCG)
702 state->flags |= (1<<estCGP);
703 if (state->cg_p == NULL)
705 /* cg_p is not stored in the tpx file, so we need to allocate it */
706 snew(state->cg_p, state->nalloc);
710 state->nnhpres = 0;
711 if (ir->ePBC != epbcNONE)
713 state->flags |= (1<<estBOX);
714 if (inputrecPreserveShape(ir))
716 state->flags |= (1<<estBOX_REL);
718 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
720 state->flags |= (1<<estBOXV);
722 if (ir->epc != epcNO)
724 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
726 state->nnhpres = 1;
727 state->flags |= (1<<estNHPRES_XI);
728 state->flags |= (1<<estNHPRES_VXI);
729 state->flags |= (1<<estSVIR_PREV);
730 state->flags |= (1<<estFVIR_PREV);
731 state->flags |= (1<<estVETA);
732 state->flags |= (1<<estVOL0);
734 else
736 state->flags |= (1<<estPRES_PREV);
741 if (ir->etc == etcNOSEHOOVER)
743 state->flags |= (1<<estNH_XI);
744 state->flags |= (1<<estNH_VXI);
747 if (ir->etc == etcVRESCALE)
749 state->flags |= (1<<estTC_INT);
752 init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
753 init_ekinstate(&state->ekinstate, ir);
754 snew(state->enerhist, 1);
755 init_energyhistory(state->enerhist);
756 init_df_history(&state->dfhist, ir->fepvals->n_lambda);
757 state->swapstate.eSwapCoords = ir->eSwapCoords;