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.
40 #include "gromacs/legacyheaders/md_support.h"
44 #include "gromacs/domdec/domdec.h"
45 #include "gromacs/fileio/trx.h"
46 #include "gromacs/legacyheaders/md_logging.h"
47 #include "gromacs/legacyheaders/mdrun.h"
48 #include "gromacs/legacyheaders/names.h"
49 #include "gromacs/legacyheaders/nrnb.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/vcm.h"
52 #include "gromacs/legacyheaders/types/commrec.h"
53 #include "gromacs/legacyheaders/types/group.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/mdlib/mdrun_signalling.h"
56 #include "gromacs/timing/wallcycle.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/utility/arrayref.h"
59 #include "gromacs/utility/cstringutil.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/snprintf.h"
64 /* check which of the multisim simulations has the shortest number of
65 steps and return that number of nsteps */
66 gmx_int64_t
get_multisim_nsteps(const t_commrec
*cr
,
69 gmx_int64_t steps_out
;
76 snew(buf
, cr
->ms
->nsim
);
78 buf
[cr
->ms
->sim
] = nsteps
;
79 gmx_sumli_sim(cr
->ms
->nsim
, buf
, cr
->ms
);
82 for (s
= 0; s
< cr
->ms
->nsim
; s
++)
84 /* find the smallest positive number */
85 if (buf
[s
] >= 0 && ((steps_out
< 0) || (buf
[s
] < steps_out
)) )
92 /* if we're the limiting simulation, don't do anything */
93 if (steps_out
>= 0 && steps_out
< nsteps
)
96 snprintf(strbuf
, 255, "Will stop simulation %%d after %s steps (another simulation will end then).\n", "%" GMX_PRId64
);
97 fprintf(stderr
, strbuf
, cr
->ms
->sim
, steps_out
);
100 /* broadcast to non-masters */
101 gmx_bcast(sizeof(gmx_int64_t
), &steps_out
, cr
);
105 int multisim_min(const gmx_multisim_t
*ms
, int nmin
, int n
)
108 gmx_bool bPos
, bEqual
;
113 gmx_sumi_sim(ms
->nsim
, buf
, ms
);
116 for (s
= 0; s
< ms
->nsim
; s
++)
118 bPos
= bPos
&& (buf
[s
] > 0);
119 bEqual
= bEqual
&& (buf
[s
] == buf
[0]);
125 nmin
= std::min(nmin
, buf
[0]);
129 /* Find the least common multiple */
130 for (d
= 2; d
< nmin
; d
++)
133 while (s
< ms
->nsim
&& d
% buf
[s
] == 0)
139 /* We found the LCM and it is less than nmin */
151 int multisim_nstsimsync(const t_commrec
*cr
,
152 const t_inputrec
*ir
, int repl_ex_nst
)
159 nmin
= multisim_min(cr
->ms
, nmin
, ir
->nstlist
);
160 nmin
= multisim_min(cr
->ms
, nmin
, ir
->nstcalcenergy
);
161 nmin
= multisim_min(cr
->ms
, nmin
, repl_ex_nst
);
164 gmx_fatal(FARGS
, "Can not find an appropriate interval for inter-simulation communication, since nstlist, nstcalcenergy and -replex are all <= 0");
166 /* Avoid inter-simulation communication at every (second) step */
173 gmx_bcast(sizeof(int), &nmin
, cr
);
178 void copy_coupling_state(t_state
*statea
, t_state
*stateb
,
179 gmx_ekindata_t
*ekinda
, gmx_ekindata_t
*ekindb
, t_grpopts
* opts
)
182 /* MRS note -- might be able to get rid of some of the arguments. Look over it when it's all debugged */
186 /* Make sure we have enough space for x and v */
187 if (statea
->nalloc
> stateb
->nalloc
)
189 stateb
->nalloc
= statea
->nalloc
;
190 srenew(stateb
->x
, stateb
->nalloc
);
191 srenew(stateb
->v
, stateb
->nalloc
);
194 stateb
->natoms
= statea
->natoms
;
195 stateb
->ngtc
= statea
->ngtc
;
196 stateb
->nnhpres
= statea
->nnhpres
;
197 stateb
->veta
= statea
->veta
;
200 copy_mat(ekinda
->ekin
, ekindb
->ekin
);
201 for (i
= 0; i
< stateb
->ngtc
; i
++)
203 ekindb
->tcstat
[i
].T
= ekinda
->tcstat
[i
].T
;
204 ekindb
->tcstat
[i
].Th
= ekinda
->tcstat
[i
].Th
;
205 copy_mat(ekinda
->tcstat
[i
].ekinh
, ekindb
->tcstat
[i
].ekinh
);
206 copy_mat(ekinda
->tcstat
[i
].ekinf
, ekindb
->tcstat
[i
].ekinf
);
207 ekindb
->tcstat
[i
].ekinscalef_nhc
= ekinda
->tcstat
[i
].ekinscalef_nhc
;
208 ekindb
->tcstat
[i
].ekinscaleh_nhc
= ekinda
->tcstat
[i
].ekinscaleh_nhc
;
209 ekindb
->tcstat
[i
].vscale_nhc
= ekinda
->tcstat
[i
].vscale_nhc
;
212 copy_rvecn(statea
->x
, stateb
->x
, 0, stateb
->natoms
);
213 copy_rvecn(statea
->v
, stateb
->v
, 0, stateb
->natoms
);
214 copy_mat(statea
->box
, stateb
->box
);
215 copy_mat(statea
->box_rel
, stateb
->box_rel
);
216 copy_mat(statea
->boxv
, stateb
->boxv
);
218 for (i
= 0; i
< stateb
->ngtc
; i
++)
220 nc
= i
*opts
->nhchainlength
;
221 for (j
= 0; j
< opts
->nhchainlength
; j
++)
223 stateb
->nosehoover_xi
[nc
+j
] = statea
->nosehoover_xi
[nc
+j
];
224 stateb
->nosehoover_vxi
[nc
+j
] = statea
->nosehoover_vxi
[nc
+j
];
227 if (stateb
->nhpres_xi
!= NULL
)
229 for (i
= 0; i
< stateb
->nnhpres
; i
++)
231 nc
= i
*opts
->nhchainlength
;
232 for (j
= 0; j
< opts
->nhchainlength
; j
++)
234 stateb
->nhpres_xi
[nc
+j
] = statea
->nhpres_xi
[nc
+j
];
235 stateb
->nhpres_vxi
[nc
+j
] = statea
->nhpres_vxi
[nc
+j
];
241 real
compute_conserved_from_auxiliary(t_inputrec
*ir
, t_state
*state
, t_extmass
*MassQ
)
251 quantity
= NPT_energy(ir
, state
, MassQ
);
254 quantity
= vrescale_energy(&(ir
->opts
), state
->therm_integral
);
262 void compute_globals(FILE *fplog
, gmx_global_stat_t gstat
, t_commrec
*cr
, t_inputrec
*ir
,
263 t_forcerec
*fr
, gmx_ekindata_t
*ekind
,
264 t_state
*state
, t_state
*state_global
, t_mdatoms
*mdatoms
,
265 t_nrnb
*nrnb
, t_vcm
*vcm
, gmx_wallcycle_t wcycle
,
266 gmx_enerdata_t
*enerd
, tensor force_vir
, tensor shake_vir
, tensor total_vir
,
267 tensor pres
, rvec mu_tot
, gmx_constr_t constr
,
268 struct gmx_signalling_t
*gs
, gmx_bool bInterSimGS
,
269 matrix box
, gmx_mtop_t
*top_global
,
270 gmx_bool
*bSumEkinhOld
, int flags
)
272 tensor corr_vir
, corr_pres
;
273 gmx_bool bEner
, bPres
, bTemp
;
274 gmx_bool bStopCM
, bGStat
,
275 bReadEkin
, bEkinAveVel
, bScaleEkin
, bConstrain
;
276 real prescorr
, enercorr
, dvdlcorr
, dvdl_ekin
;
278 /* translate CGLO flags to gmx_booleans */
279 bStopCM
= flags
& CGLO_STOPCM
;
280 bGStat
= flags
& CGLO_GSTAT
;
281 bReadEkin
= (flags
& CGLO_READEKIN
);
282 bScaleEkin
= (flags
& CGLO_SCALEEKIN
);
283 bEner
= flags
& CGLO_ENERGY
;
284 bTemp
= flags
& CGLO_TEMPERATURE
;
285 bPres
= (flags
& CGLO_PRESSURE
);
286 bConstrain
= (flags
& CGLO_CONSTRAINT
);
288 /* we calculate a full state kinetic energy either with full-step velocity verlet
289 or half step where we need the pressure */
291 bEkinAveVel
= (ir
->eI
== eiVV
|| (ir
->eI
== eiVVAK
&& bPres
) || bReadEkin
);
293 /* in initalization, it sums the shake virial in vv, and to
294 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
296 /* ########## Kinetic energy ############## */
300 /* Non-equilibrium MD: this is parallellized, but only does communication
301 * when there really is NEMD.
304 if (PAR(cr
) && (ekind
->bNEMD
))
306 accumulate_u(cr
, &(ir
->opts
), ekind
);
311 restore_ekinstate_from_state(cr
, ekind
, &state_global
->ekinstate
);
315 calc_ke_part(state
, &(ir
->opts
), mdatoms
, ekind
, nrnb
, bEkinAveVel
);
321 /* Calculate center of mass velocity if necessary, also parallellized */
324 calc_vcm_grp(0, mdatoms
->homenr
, mdatoms
,
325 state
->x
, state
->v
, vcm
);
328 if (bTemp
|| bStopCM
|| bPres
|| bEner
|| bConstrain
)
332 /* We will not sum ekinh_old,
333 * so signal that we still have to do it.
335 *bSumEkinhOld
= TRUE
;
340 gmx::ArrayRef
<real
> signalBuffer
= prepareSignalBuffer(gs
);
343 wallcycle_start(wcycle
, ewcMoveE
);
344 global_stat(fplog
, gstat
, cr
, enerd
, force_vir
, shake_vir
, mu_tot
,
345 ir
, ekind
, constr
, bStopCM
? vcm
: NULL
,
346 signalBuffer
.size(), signalBuffer
.data(),
348 *bSumEkinhOld
, flags
);
349 wallcycle_stop(wcycle
, ewcMoveE
);
351 handleSignals(gs
, cr
, bInterSimGS
);
352 *bSumEkinhOld
= FALSE
;
356 if (!ekind
->bNEMD
&& debug
&& bTemp
&& (vcm
->nr
> 0))
360 state
->v
, vcm
->group_p
[0],
361 mdatoms
->massT
, mdatoms
->tmass
, ekind
->ekin
);
364 /* Do center of mass motion removal */
367 check_cm_grp(fplog
, vcm
, ir
, 1);
368 do_stopcm_grp(0, mdatoms
->homenr
, mdatoms
->cVCM
,
369 state
->x
, state
->v
, vcm
);
370 inc_nrnb(nrnb
, eNR_STOPCM
, mdatoms
->homenr
);
375 /* Calculate the amplitude of the cosine velocity profile */
376 ekind
->cosacc
.vcos
= ekind
->cosacc
.mvcos
/mdatoms
->tmass
;
381 /* Sum the kinetic energies of the groups & calc temp */
382 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with IR_NPT_TROTTER */
383 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
384 Leap with AveVel is not supported; it's not clear that it will actually work.
385 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
386 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
388 enerd
->term
[F_TEMP
] = sum_ekin(&(ir
->opts
), ekind
, &dvdl_ekin
,
389 bEkinAveVel
, bScaleEkin
);
390 enerd
->dvdl_lin
[efptMASS
] = (double) dvdl_ekin
;
392 enerd
->term
[F_EKIN
] = trace(ekind
->ekin
);
395 /* ########## Long range energy information ###### */
397 if (bEner
|| bPres
|| bConstrain
)
399 calc_dispcorr(ir
, fr
, top_global
->natoms
, box
, state
->lambda
[efptVDW
],
400 corr_pres
, corr_vir
, &prescorr
, &enercorr
, &dvdlcorr
);
405 enerd
->term
[F_DISPCORR
] = enercorr
;
406 enerd
->term
[F_EPOT
] += enercorr
;
407 enerd
->term
[F_DVDL_VDW
] += dvdlcorr
;
410 /* ########## Now pressure ############## */
411 if (bPres
|| bConstrain
)
414 m_add(force_vir
, shake_vir
, total_vir
);
416 /* Calculate pressure and apply LR correction if PPPM is used.
417 * Use the box from last timestep since we already called update().
420 enerd
->term
[F_PRES
] = calc_pres(fr
->ePBC
, ir
->nwall
, box
, ekind
->ekin
, total_vir
, pres
);
422 /* Calculate long range corrections to pressure and energy */
423 /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
424 and computes enerd->term[F_DISPCORR]. Also modifies the
425 total_vir and pres tesors */
427 m_add(total_vir
, corr_vir
, total_vir
);
428 m_add(pres
, corr_pres
, pres
);
429 enerd
->term
[F_PDISPCORR
] = prescorr
;
430 enerd
->term
[F_PRES
] += prescorr
;
434 void check_nst_param(FILE *fplog
, t_commrec
*cr
,
435 const char *desc_nst
, int nst
,
436 const char *desc_p
, int *p
)
438 if (*p
> 0 && *p
% nst
!= 0)
440 /* Round up to the next multiple of nst */
441 *p
= ((*p
)/nst
+ 1)*nst
;
442 md_print_warn(cr
, fplog
,
443 "NOTE: %s changes %s to %d\n", desc_nst
, desc_p
, *p
);
447 void set_current_lambdas(gmx_int64_t step
, t_lambda
*fepvals
, gmx_bool bRerunMD
,
448 t_trxframe
*rerun_fr
, t_state
*state_global
, t_state
*state
, double lam0
[])
449 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
450 requiring different logic. */
453 int i
, fep_state
= 0;
456 if (rerun_fr
->bLambda
)
458 if (fepvals
->delta_lambda
== 0)
460 state_global
->lambda
[efptFEP
] = rerun_fr
->lambda
;
461 for (i
= 0; i
< efptNR
; i
++)
465 state
->lambda
[i
] = state_global
->lambda
[i
];
471 /* find out between which two value of lambda we should be */
472 frac
= (step
*fepvals
->delta_lambda
);
473 fep_state
= static_cast<int>(floor(frac
*fepvals
->n_lambda
));
474 /* interpolate between this state and the next */
475 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
476 frac
= (frac
*fepvals
->n_lambda
)-fep_state
;
477 for (i
= 0; i
< efptNR
; i
++)
479 state_global
->lambda
[i
] = lam0
[i
] + (fepvals
->all_lambda
[i
][fep_state
]) +
480 frac
*(fepvals
->all_lambda
[i
][fep_state
+1]-fepvals
->all_lambda
[i
][fep_state
]);
484 else if (rerun_fr
->bFepState
)
486 state_global
->fep_state
= rerun_fr
->fep_state
;
487 for (i
= 0; i
< efptNR
; i
++)
489 state_global
->lambda
[i
] = fepvals
->all_lambda
[i
][fep_state
];
495 if (fepvals
->delta_lambda
!= 0)
497 /* find out between which two value of lambda we should be */
498 frac
= (step
*fepvals
->delta_lambda
);
499 if (fepvals
->n_lambda
> 0)
501 fep_state
= static_cast<int>(floor(frac
*fepvals
->n_lambda
));
502 /* interpolate between this state and the next */
503 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
504 frac
= (frac
*fepvals
->n_lambda
)-fep_state
;
505 for (i
= 0; i
< efptNR
; i
++)
507 state_global
->lambda
[i
] = lam0
[i
] + (fepvals
->all_lambda
[i
][fep_state
]) +
508 frac
*(fepvals
->all_lambda
[i
][fep_state
+1]-fepvals
->all_lambda
[i
][fep_state
]);
513 for (i
= 0; i
< efptNR
; i
++)
515 state_global
->lambda
[i
] = lam0
[i
] + frac
;
521 if (state
->fep_state
> 0)
523 state_global
->fep_state
= state
->fep_state
; /* state->fep is the one updated by bExpanded */
524 for (i
= 0; i
< efptNR
; i
++)
526 state_global
->lambda
[i
] = fepvals
->all_lambda
[i
][state_global
->fep_state
];
531 for (i
= 0; i
< efptNR
; i
++)
533 state
->lambda
[i
] = state_global
->lambda
[i
];
537 static void min_zero(int *n
, int i
)
539 if (i
> 0 && (*n
== 0 || i
< *n
))
545 static int lcd4(int i1
, int i2
, int i3
, int i4
)
556 gmx_incons("All 4 inputs for determininig nstglobalcomm are <= 0");
559 while (nst
> 1 && ((i1
> 0 && i1
% nst
!= 0) ||
560 (i2
> 0 && i2
% nst
!= 0) ||
561 (i3
> 0 && i3
% nst
!= 0) ||
562 (i4
> 0 && i4
% nst
!= 0)))
570 int check_nstglobalcomm(FILE *fplog
, t_commrec
*cr
,
571 int nstglobalcomm
, t_inputrec
*ir
)
573 if (!EI_DYNAMICS(ir
->eI
))
578 if (nstglobalcomm
== -1)
580 if (!(ir
->nstcalcenergy
> 0 ||
586 if (ir
->nstenergy
> 0 && ir
->nstenergy
< nstglobalcomm
)
588 nstglobalcomm
= ir
->nstenergy
;
593 /* Ensure that we do timely global communication for
594 * (possibly) each of the four following options.
596 nstglobalcomm
= lcd4(ir
->nstcalcenergy
,
598 ir
->etc
!= etcNO
? ir
->nsttcouple
: 0,
599 ir
->epc
!= epcNO
? ir
->nstpcouple
: 0);
604 if (ir
->nstlist
> 0 &&
605 nstglobalcomm
> ir
->nstlist
&& nstglobalcomm
% ir
->nstlist
!= 0)
607 nstglobalcomm
= (nstglobalcomm
/ ir
->nstlist
)*ir
->nstlist
;
608 md_print_warn(cr
, fplog
, "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d\n", nstglobalcomm
);
610 if (ir
->nstcalcenergy
> 0)
612 check_nst_param(fplog
, cr
, "-gcom", nstglobalcomm
,
613 "nstcalcenergy", &ir
->nstcalcenergy
);
615 if (ir
->etc
!= etcNO
&& ir
->nsttcouple
> 0)
617 check_nst_param(fplog
, cr
, "-gcom", nstglobalcomm
,
618 "nsttcouple", &ir
->nsttcouple
);
620 if (ir
->epc
!= epcNO
&& ir
->nstpcouple
> 0)
622 check_nst_param(fplog
, cr
, "-gcom", nstglobalcomm
,
623 "nstpcouple", &ir
->nstpcouple
);
626 check_nst_param(fplog
, cr
, "-gcom", nstglobalcomm
,
627 "nstenergy", &ir
->nstenergy
);
629 check_nst_param(fplog
, cr
, "-gcom", nstglobalcomm
,
630 "nstlog", &ir
->nstlog
);
633 if (ir
->comm_mode
!= ecmNO
&& ir
->nstcomm
< nstglobalcomm
)
635 md_print_warn(cr
, fplog
, "WARNING: Changing nstcomm from %d to %d\n",
636 ir
->nstcomm
, nstglobalcomm
);
637 ir
->nstcomm
= nstglobalcomm
;
640 return nstglobalcomm
;
643 void check_ir_old_tpx_versions(t_commrec
*cr
, FILE *fplog
,
644 t_inputrec
*ir
, gmx_mtop_t
*mtop
)
646 /* Check required for old tpx files */
647 if (IR_TWINRANGE(*ir
) && ir
->nstlist
> 1 &&
648 ir
->nstcalcenergy
% ir
->nstlist
!= 0)
650 md_print_warn(cr
, fplog
, "Old tpr file with twin-range settings: modifying energy calculation and/or T/P-coupling frequencies\n");
652 if (gmx_mtop_ftype_count(mtop
, F_CONSTR
) +
653 gmx_mtop_ftype_count(mtop
, F_CONSTRNC
) > 0 &&
654 ir
->eConstrAlg
== econtSHAKE
)
656 md_print_warn(cr
, fplog
, "With twin-range cut-off's and SHAKE the virial and pressure are incorrect\n");
657 if (ir
->epc
!= epcNO
)
659 gmx_fatal(FARGS
, "Can not do pressure coupling with twin-range cut-off's and SHAKE");
662 check_nst_param(fplog
, cr
, "nstlist", ir
->nstlist
,
663 "nstcalcenergy", &ir
->nstcalcenergy
);
664 if (ir
->epc
!= epcNO
)
666 check_nst_param(fplog
, cr
, "nstlist", ir
->nstlist
,
667 "nstpcouple", &ir
->nstpcouple
);
669 check_nst_param(fplog
, cr
, "nstcalcenergy", ir
->nstcalcenergy
,
670 "nstenergy", &ir
->nstenergy
);
671 check_nst_param(fplog
, cr
, "nstcalcenergy", ir
->nstcalcenergy
,
672 "nstlog", &ir
->nstlog
);
673 if (ir
->efep
!= efepNO
)
675 check_nst_param(fplog
, cr
, "nstcalcenergy", ir
->nstcalcenergy
,
676 "nstdhdl", &ir
->fepvals
->nstdhdl
);
680 if (EI_VV(ir
->eI
) && IR_TWINRANGE(*ir
) && ir
->nstlist
> 1)
682 gmx_fatal(FARGS
, "Twin-range multiple time stepping does not work with integrator %s.", ei_names
[ir
->eI
]);
686 void rerun_parallel_comm(t_commrec
*cr
, t_trxframe
*fr
,
687 gmx_bool
*bNotLastFrame
)
691 if (MASTER(cr
) && !*bNotLastFrame
)
697 gmx_bcast(sizeof(*fr
), fr
, cr
);
701 *bNotLastFrame
= (fr
->natoms
>= 0);
705 void set_state_entries(t_state
*state
, const t_inputrec
*ir
)
707 /* The entries in the state in the tpx file might not correspond
708 * with what is needed, so we correct this here.
711 if (ir
->efep
!= efepNO
|| ir
->bExpanded
)
713 state
->flags
|= (1<<estLAMBDA
);
714 state
->flags
|= (1<<estFEPSTATE
);
716 state
->flags
|= (1<<estX
);
717 if (state
->lambda
== NULL
)
719 snew(state
->lambda
, efptNR
);
721 if (state
->x
== NULL
)
723 snew(state
->x
, state
->nalloc
);
725 if (EI_DYNAMICS(ir
->eI
))
727 state
->flags
|= (1<<estV
);
728 if (state
->v
== NULL
)
730 snew(state
->v
, state
->nalloc
);
735 state
->flags
|= (1<<estSDX
);
736 if (state
->sd_X
== NULL
)
738 /* sd_X is not stored in the tpx file, so we need to allocate it */
739 snew(state
->sd_X
, state
->nalloc
);
744 state
->flags
|= (1<<estCGP
);
745 if (state
->cg_p
== NULL
)
747 /* cg_p is not stored in the tpx file, so we need to allocate it */
748 snew(state
->cg_p
, state
->nalloc
);
753 if (ir
->ePBC
!= epbcNONE
)
755 state
->flags
|= (1<<estBOX
);
756 if (PRESERVE_SHAPE(*ir
))
758 state
->flags
|= (1<<estBOX_REL
);
760 if ((ir
->epc
== epcPARRINELLORAHMAN
) || (ir
->epc
== epcMTTK
))
762 state
->flags
|= (1<<estBOXV
);
764 if (ir
->epc
!= epcNO
)
766 if (IR_NPT_TROTTER(ir
) || (IR_NPH_TROTTER(ir
)))
769 state
->flags
|= (1<<estNHPRES_XI
);
770 state
->flags
|= (1<<estNHPRES_VXI
);
771 state
->flags
|= (1<<estSVIR_PREV
);
772 state
->flags
|= (1<<estFVIR_PREV
);
773 state
->flags
|= (1<<estVETA
);
774 state
->flags
|= (1<<estVOL0
);
778 state
->flags
|= (1<<estPRES_PREV
);
783 if (ir
->etc
== etcNOSEHOOVER
)
785 state
->flags
|= (1<<estNH_XI
);
786 state
->flags
|= (1<<estNH_VXI
);
789 if (ir
->etc
== etcVRESCALE
)
791 state
->flags
|= (1<<estTC_INT
);
794 init_gtc_state(state
, state
->ngtc
, state
->nnhpres
, ir
->opts
.nhchainlength
); /* allocate the space for nose-hoover chains */
795 init_ekinstate(&state
->ekinstate
, ir
);
797 init_energyhistory(&state
->enerhist
);
798 init_df_history(&state
->dfhist
, ir
->fepvals
->n_lambda
);
799 state
->swapstate
.eSwapCoords
= ir
->eSwapCoords
;