Introduce ObservablesHistory container
[gromacs.git] / src / gromacs / mdlib / md_support.cpp
blob61224b1bc984d8a92b09bf5faeddb7ec5b898c6e
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,2016,2017, 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 <climits>
44 #include <algorithm>
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/gmxlib/nrnb.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/mdrun.h"
51 #include "gromacs/mdlib/sim_util.h"
52 #include "gromacs/mdlib/simulationsignal.h"
53 #include "gromacs/mdlib/tgroup.h"
54 #include "gromacs/mdlib/update.h"
55 #include "gromacs/mdlib/vcm.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/mdtypes/df_history.h"
58 #include "gromacs/mdtypes/energyhistory.h"
59 #include "gromacs/mdtypes/group.h"
60 #include "gromacs/mdtypes/inputrec.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/mdtypes/state.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/timing/wallcycle.h"
65 #include "gromacs/topology/mtop_util.h"
66 #include "gromacs/trajectory/trajectoryframe.h"
67 #include "gromacs/utility/arrayref.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/gmxassert.h"
71 #include "gromacs/utility/logger.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/snprintf.h"
75 // TODO move this to multi-sim module
76 bool multisim_int_all_are_equal(const gmx_multisim_t *ms,
77 gmx_int64_t value)
79 bool allValuesAreEqual = true;
80 gmx_int64_t *buf;
82 GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
84 snew(buf, ms->nsim);
85 /* send our value to all other master ranks, receive all of theirs */
86 buf[ms->sim] = value;
87 gmx_sumli_sim(ms->nsim, buf, ms);
89 for (int s = 0; s < ms->nsim; s++)
91 if (buf[s] != value)
93 allValuesAreEqual = false;
94 break;
98 sfree(buf);
100 return allValuesAreEqual;
103 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
105 int *buf;
106 gmx_bool bPos, bEqual;
107 int s, d;
109 snew(buf, ms->nsim);
110 buf[ms->sim] = n;
111 gmx_sumi_sim(ms->nsim, buf, ms);
112 bPos = TRUE;
113 bEqual = TRUE;
114 for (s = 0; s < ms->nsim; s++)
116 bPos = bPos && (buf[s] > 0);
117 bEqual = bEqual && (buf[s] == buf[0]);
119 if (bPos)
121 if (bEqual)
123 nmin = std::min(nmin, buf[0]);
125 else
127 /* Find the least common multiple */
128 for (d = 2; d < nmin; d++)
130 s = 0;
131 while (s < ms->nsim && d % buf[s] == 0)
133 s++;
135 if (s == ms->nsim)
137 /* We found the LCM and it is less than nmin */
138 nmin = d;
139 break;
144 sfree(buf);
146 return nmin;
149 /* TODO Specialize this routine into init-time and loop-time versions?
150 e.g. bReadEkin is only true when restoring from checkpoint */
151 void compute_globals(FILE *fplog, gmx_global_stat *gstat, t_commrec *cr, t_inputrec *ir,
152 t_forcerec *fr, gmx_ekindata_t *ekind,
153 t_state *state, t_mdatoms *mdatoms,
154 t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
155 gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
156 tensor pres, rvec mu_tot, gmx_constr_t constr,
157 gmx::SimulationSignaller *signalCoordinator,
158 matrix box, int *totalNumberOfBondedInteractions,
159 gmx_bool *bSumEkinhOld, int flags)
161 tensor corr_vir, corr_pres;
162 gmx_bool bEner, bPres, bTemp;
163 gmx_bool bStopCM, bGStat,
164 bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
165 real prescorr, enercorr, dvdlcorr, dvdl_ekin;
167 /* translate CGLO flags to gmx_booleans */
168 bStopCM = flags & CGLO_STOPCM;
169 bGStat = flags & CGLO_GSTAT;
170 bReadEkin = (flags & CGLO_READEKIN);
171 bScaleEkin = (flags & CGLO_SCALEEKIN);
172 bEner = flags & CGLO_ENERGY;
173 bTemp = flags & CGLO_TEMPERATURE;
174 bPres = (flags & CGLO_PRESSURE);
175 bConstrain = (flags & CGLO_CONSTRAINT);
177 /* we calculate a full state kinetic energy either with full-step velocity verlet
178 or half step where we need the pressure */
180 bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
182 /* in initalization, it sums the shake virial in vv, and to
183 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
185 /* ########## Kinetic energy ############## */
187 if (bTemp)
189 /* Non-equilibrium MD: this is parallellized, but only does communication
190 * when there really is NEMD.
193 if (PAR(cr) && (ekind->bNEMD))
195 accumulate_u(cr, &(ir->opts), ekind);
197 if (!bReadEkin)
199 calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
203 /* Calculate center of mass velocity if necessary, also parallellized */
204 if (bStopCM)
206 calc_vcm_grp(0, mdatoms->homenr, mdatoms,
207 as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), vcm);
210 if (bTemp || bStopCM || bPres || bEner || bConstrain)
212 if (!bGStat)
214 /* We will not sum ekinh_old,
215 * so signal that we still have to do it.
217 *bSumEkinhOld = TRUE;
220 else
222 gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
223 if (PAR(cr))
225 wallcycle_start(wcycle, ewcMoveE);
226 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
227 ir, ekind, constr, bStopCM ? vcm : nullptr,
228 signalBuffer.size(), signalBuffer.data(),
229 totalNumberOfBondedInteractions,
230 *bSumEkinhOld, flags);
231 wallcycle_stop(wcycle, ewcMoveE);
233 signalCoordinator->finalizeSignals();
234 *bSumEkinhOld = FALSE;
238 if (!ekind->bNEMD && debug && bTemp && (vcm->nr > 0))
240 correct_ekin(debug,
241 0, mdatoms->homenr,
242 as_rvec_array(state->v.data()), vcm->group_p[0],
243 mdatoms->massT, mdatoms->tmass, ekind->ekin);
246 /* Do center of mass motion removal */
247 if (bStopCM)
249 check_cm_grp(fplog, vcm, ir, 1);
250 do_stopcm_grp(0, mdatoms->homenr, mdatoms->cVCM,
251 as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), vcm);
252 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
255 if (bEner)
257 /* Calculate the amplitude of the cosine velocity profile */
258 ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
261 if (bTemp)
263 /* Sum the kinetic energies of the groups & calc temp */
264 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
265 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
266 Leap with AveVel is not supported; it's not clear that it will actually work.
267 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
268 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
270 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
271 bEkinAveVel, bScaleEkin);
272 enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
274 enerd->term[F_EKIN] = trace(ekind->ekin);
277 /* ########## Long range energy information ###### */
279 if (bEner || bPres || bConstrain)
281 calc_dispcorr(ir, fr, box, state->lambda[efptVDW],
282 corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
285 if (bEner)
287 enerd->term[F_DISPCORR] = enercorr;
288 enerd->term[F_EPOT] += enercorr;
289 enerd->term[F_DVDL_VDW] += dvdlcorr;
292 /* ########## Now pressure ############## */
293 if (bPres || bConstrain)
296 m_add(force_vir, shake_vir, total_vir);
298 /* Calculate pressure and apply LR correction if PPPM is used.
299 * Use the box from last timestep since we already called update().
302 enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
304 /* Calculate long range corrections to pressure and energy */
305 /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
306 and computes enerd->term[F_DISPCORR]. Also modifies the
307 total_vir and pres tesors */
309 m_add(total_vir, corr_vir, total_vir);
310 m_add(pres, corr_pres, pres);
311 enerd->term[F_PDISPCORR] = prescorr;
312 enerd->term[F_PRES] += prescorr;
316 /* check whether an 'nst'-style parameter p is a multiple of nst, and
317 set it to be one if not, with a warning. */
318 static void check_nst_param(const gmx::MDLogger &mdlog,
319 const char *desc_nst, int nst,
320 const char *desc_p, int *p)
322 if (*p > 0 && *p % nst != 0)
324 /* Round up to the next multiple of nst */
325 *p = ((*p)/nst + 1)*nst;
326 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
327 "NOTE: %s changes %s to %d", desc_nst, desc_p, *p);
331 void set_current_lambdas(gmx_int64_t step, t_lambda *fepvals, gmx_bool bRerunMD,
332 t_trxframe *rerun_fr, t_state *state_global, t_state *state, double lam0[])
333 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
334 requiring different logic. */
336 real frac;
337 int i, fep_state = 0;
338 if (bRerunMD)
340 if (rerun_fr->bLambda)
342 if (fepvals->delta_lambda == 0)
344 state_global->lambda[efptFEP] = rerun_fr->lambda;
345 for (i = 0; i < efptNR; i++)
347 if (i != efptFEP)
349 state->lambda[i] = state_global->lambda[i];
353 else
355 /* find out between which two value of lambda we should be */
356 frac = (step*fepvals->delta_lambda);
357 fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
358 /* interpolate between this state and the next */
359 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
360 frac = (frac*fepvals->n_lambda)-fep_state;
361 for (i = 0; i < efptNR; i++)
363 state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
364 frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
368 else if (rerun_fr->bFepState)
370 state_global->fep_state = rerun_fr->fep_state;
371 for (i = 0; i < efptNR; i++)
373 state_global->lambda[i] = fepvals->all_lambda[i][fep_state];
377 else
379 if (fepvals->delta_lambda != 0)
381 /* find out between which two value of lambda we should be */
382 frac = (step*fepvals->delta_lambda);
383 if (fepvals->n_lambda > 0)
385 fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
386 /* interpolate between this state and the next */
387 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
388 frac = (frac*fepvals->n_lambda)-fep_state;
389 for (i = 0; i < efptNR; i++)
391 state_global->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
392 frac*(fepvals->all_lambda[i][fep_state+1]-fepvals->all_lambda[i][fep_state]);
395 else
397 for (i = 0; i < efptNR; i++)
399 state_global->lambda[i] = lam0[i] + frac;
403 else
405 /* if < 0, fep_state was never defined, and we should not set lambda from the state */
406 if (state_global->fep_state > -1)
408 state_global->fep_state = state->fep_state; /* state->fep_state is the one updated by bExpanded */
409 for (i = 0; i < efptNR; i++)
411 state_global->lambda[i] = fepvals->all_lambda[i][state_global->fep_state];
416 for (i = 0; i < efptNR; i++)
418 state->lambda[i] = state_global->lambda[i];
422 static void min_zero(int *n, int i)
424 if (i > 0 && (*n == 0 || i < *n))
426 *n = i;
430 static int lcd4(int i1, int i2, int i3, int i4)
432 int nst;
434 nst = 0;
435 min_zero(&nst, i1);
436 min_zero(&nst, i2);
437 min_zero(&nst, i3);
438 min_zero(&nst, i4);
439 if (nst == 0)
441 gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
444 while (nst > 1 && ((i1 > 0 && i1 % nst != 0) ||
445 (i2 > 0 && i2 % nst != 0) ||
446 (i3 > 0 && i3 % nst != 0) ||
447 (i4 > 0 && i4 % nst != 0)))
449 nst--;
452 return nst;
455 int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputrec *ir)
457 if (!EI_DYNAMICS(ir->eI))
459 nstglobalcomm = 1;
462 if (nstglobalcomm == -1)
464 // Set up the default behaviour
465 if (!(ir->nstcalcenergy > 0 ||
466 ir->nstlist > 0 ||
467 ir->etc != etcNO ||
468 ir->epc != epcNO))
470 /* The user didn't choose the period for anything
471 important, so we just make sure we can send signals and
472 write output suitably. */
473 nstglobalcomm = 10;
474 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
476 nstglobalcomm = ir->nstenergy;
479 else
481 /* The user has made a choice (perhaps implicitly), so we
482 * ensure that we do timely intra-simulation communication
483 * for (possibly) each of the four parts that care.
485 * TODO Does the Verlet scheme (+ DD) need any
486 * communication at nstlist steps? Is the use of nstlist
487 * here a leftover of the twin-range scheme? Can we remove
488 * nstlist when we remove the group scheme?
490 nstglobalcomm = lcd4(ir->nstcalcenergy,
491 ir->nstlist,
492 ir->etc != etcNO ? ir->nsttcouple : 0,
493 ir->epc != epcNO ? ir->nstpcouple : 0);
496 else
498 // Check that the user's choice of mdrun -gcom will work
499 if (ir->nstlist > 0 &&
500 nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
502 nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
503 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
504 "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d",
505 nstglobalcomm);
507 if (ir->nstcalcenergy > 0)
509 check_nst_param(mdlog, "-gcom", nstglobalcomm,
510 "nstcalcenergy", &ir->nstcalcenergy);
512 if (ir->etc != etcNO && ir->nsttcouple > 0)
514 check_nst_param(mdlog, "-gcom", nstglobalcomm,
515 "nsttcouple", &ir->nsttcouple);
517 if (ir->epc != epcNO && ir->nstpcouple > 0)
519 check_nst_param(mdlog, "-gcom", nstglobalcomm,
520 "nstpcouple", &ir->nstpcouple);
523 check_nst_param(mdlog, "-gcom", nstglobalcomm,
524 "nstenergy", &ir->nstenergy);
526 check_nst_param(mdlog, "-gcom", nstglobalcomm,
527 "nstlog", &ir->nstlog);
530 if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
532 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
533 "WARNING: Changing nstcomm from %d to %d",
534 ir->nstcomm, nstglobalcomm);
535 ir->nstcomm = nstglobalcomm;
538 GMX_LOG(mdlog.info).appendTextFormatted(
539 "Intra-simulation communication will occur every %d steps.\n", nstglobalcomm);
540 return nstglobalcomm;
543 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
544 gmx_bool *bLastStep)
546 rvec *xp, *vp;
548 if (MASTER(cr) && *bLastStep)
550 fr->natoms = -1;
552 xp = fr->x;
553 vp = fr->v;
554 gmx_bcast(sizeof(*fr), fr, cr);
555 fr->x = xp;
556 fr->v = vp;
558 *bLastStep = (fr->natoms < 0);
562 // TODO Most of this logic seems to belong in the respective modules
563 void set_state_entries(t_state *state, const t_inputrec *ir)
565 /* The entries in the state in the tpx file might not correspond
566 * with what is needed, so we correct this here.
568 state->flags = 0;
569 if (ir->efep != efepNO || ir->bExpanded)
571 state->flags |= (1<<estLAMBDA);
572 state->flags |= (1<<estFEPSTATE);
574 state->flags |= (1<<estX);
575 GMX_RELEASE_ASSERT(state->x.size() >= static_cast<unsigned int>(state->natoms), "We should start a run with an initialized state->x");
576 if (EI_DYNAMICS(ir->eI))
578 state->flags |= (1<<estV);
581 state->nnhpres = 0;
582 if (ir->ePBC != epbcNONE)
584 state->flags |= (1<<estBOX);
585 if (inputrecPreserveShape(ir))
587 state->flags |= (1<<estBOX_REL);
589 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
591 state->flags |= (1<<estBOXV);
592 state->flags |= (1<<estPRES_PREV);
594 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
596 state->nnhpres = 1;
597 state->flags |= (1<<estNHPRES_XI);
598 state->flags |= (1<<estNHPRES_VXI);
599 state->flags |= (1<<estSVIR_PREV);
600 state->flags |= (1<<estFVIR_PREV);
601 state->flags |= (1<<estVETA);
602 state->flags |= (1<<estVOL0);
606 if (ir->etc == etcNOSEHOOVER)
608 state->flags |= (1<<estNH_XI);
609 state->flags |= (1<<estNH_VXI);
612 if (ir->etc == etcVRESCALE)
614 state->flags |= (1<<estTC_INT);
617 init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
618 init_ekinstate(&state->ekinstate, ir);
620 if (ir->bExpanded)
622 snew(state->dfhist, 1);
623 init_df_history(state->dfhist, ir->fepvals->n_lambda);