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,2018, 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 "md_support.h"
47 #include "gromacs/domdec/domdec.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/mdrun.h"
52 #include "gromacs/mdlib/sim_util.h"
53 #include "gromacs/mdlib/simulationsignal.h"
54 #include "gromacs/mdlib/tgroup.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/df_history.h"
59 #include "gromacs/mdtypes/energyhistory.h"
60 #include "gromacs/mdtypes/forcerec.h"
61 #include "gromacs/mdtypes/group.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/mdtypes/state.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/topology/mtop_util.h"
68 #include "gromacs/trajectory/trajectoryframe.h"
69 #include "gromacs/utility/arrayref.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/gmxassert.h"
73 #include "gromacs/utility/logger.h"
74 #include "gromacs/utility/smalloc.h"
75 #include "gromacs/utility/snprintf.h"
77 // TODO move this to multi-sim module
78 bool multisim_int_all_are_equal(const gmx_multisim_t *ms,
81 bool allValuesAreEqual = true;
84 GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
87 /* send our value to all other master ranks, receive all of theirs */
89 gmx_sumli_sim(ms->nsim, buf, ms);
91 for (int s = 0; s < ms->nsim; s++)
95 allValuesAreEqual = false;
102 return allValuesAreEqual;
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 /* TODO Specialize this routine into init-time and loop-time versions?
152 e.g. bReadEkin is only true when restoring from checkpoint */
153 void compute_globals(FILE *fplog, gmx_global_stat *gstat, t_commrec *cr, t_inputrec *ir,
154 t_forcerec *fr, gmx_ekindata_t *ekind,
155 t_state *state, t_mdatoms *mdatoms,
156 t_nrnb *nrnb, t_vcm *vcm, gmx_wallcycle_t wcycle,
157 gmx_enerdata_t *enerd, tensor force_vir, tensor shake_vir, tensor total_vir,
158 tensor pres, rvec mu_tot, gmx::Constraints *constr,
159 gmx::SimulationSignaller *signalCoordinator,
160 matrix box, int *totalNumberOfBondedInteractions,
161 gmx_bool *bSumEkinhOld, int flags)
163 tensor corr_vir, corr_pres;
164 gmx_bool bEner, bPres, bTemp;
165 gmx_bool bStopCM, bGStat,
166 bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
167 real prescorr, enercorr, dvdlcorr, dvdl_ekin;
169 /* translate CGLO flags to gmx_booleans */
170 bStopCM = flags & CGLO_STOPCM;
171 bGStat = flags & CGLO_GSTAT;
172 bReadEkin = (flags & CGLO_READEKIN);
173 bScaleEkin = (flags & CGLO_SCALEEKIN);
174 bEner = flags & CGLO_ENERGY;
175 bTemp = flags & CGLO_TEMPERATURE;
176 bPres = (flags & CGLO_PRESSURE);
177 bConstrain = (flags & CGLO_CONSTRAINT);
179 /* we calculate a full state kinetic energy either with full-step velocity verlet
180 or half step where we need the pressure */
182 bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
184 /* in initalization, it sums the shake virial in vv, and to
185 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
187 /* ########## Kinetic energy ############## */
191 /* Non-equilibrium MD: this is parallellized, but only does communication
192 * when there really is NEMD.
195 if (PAR(cr) && (ekind->bNEMD))
197 accumulate_u(cr, &(ir->opts), ekind);
201 calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
205 /* Calculate center of mass velocity if necessary, also parallellized */
208 calc_vcm_grp(0, mdatoms->homenr, mdatoms,
209 as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), vcm);
212 if (bTemp || bStopCM || bPres || bEner || bConstrain)
216 /* We will not sum ekinh_old,
217 * so signal that we still have to do it.
219 *bSumEkinhOld = TRUE;
224 gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
227 wallcycle_start(wcycle, ewcMoveE);
228 global_stat(gstat, cr, enerd, force_vir, shake_vir, mu_tot,
229 ir, ekind, constr, bStopCM ? vcm : nullptr,
230 signalBuffer.size(), signalBuffer.data(),
231 totalNumberOfBondedInteractions,
232 *bSumEkinhOld, flags);
233 wallcycle_stop(wcycle, ewcMoveE);
235 signalCoordinator->finalizeSignals();
236 *bSumEkinhOld = FALSE;
240 /* Do center of mass motion removal */
243 check_cm_grp(fplog, vcm, ir, 1);
244 /* At initialization, do not pass x with acceleration-correction mode
245 * to avoid (incorrect) correction of the initial coordinates.
247 rvec *xPtr = nullptr;
248 if (vcm->mode == ecmANGULAR || (vcm->mode == ecmLINEAR_ACCELERATION_CORRECTION && !(flags & CGLO_INITIALIZATION)))
250 xPtr = as_rvec_array(state->x.data());
252 do_stopcm_grp(mdatoms->homenr, mdatoms->cVCM,
253 xPtr, as_rvec_array(state->v.data()), *vcm);
254 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
259 /* Calculate the amplitude of the cosine velocity profile */
260 ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
265 /* Sum the kinetic energies of the groups & calc temp */
266 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
267 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
268 Leap with AveVel is not supported; it's not clear that it will actually work.
269 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
270 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
272 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
273 bEkinAveVel, bScaleEkin);
274 enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
276 enerd->term[F_EKIN] = trace(ekind->ekin);
279 /* ########## Long range energy information ###### */
281 if (bEner || bPres || bConstrain)
283 calc_dispcorr(ir, fr, box, state->lambda[efptVDW],
284 corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
289 enerd->term[F_DISPCORR] = enercorr;
290 enerd->term[F_EPOT] += enercorr;
291 enerd->term[F_DVDL_VDW] += dvdlcorr;
294 /* ########## Now pressure ############## */
295 if (bPres || bConstrain)
298 m_add(force_vir, shake_vir, total_vir);
300 /* Calculate pressure and apply LR correction if PPPM is used.
301 * Use the box from last timestep since we already called update().
304 enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
306 /* Calculate long range corrections to pressure and energy */
307 /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
308 and computes enerd->term[F_DISPCORR]. Also modifies the
309 total_vir and pres tesors */
311 m_add(total_vir, corr_vir, total_vir);
312 m_add(pres, corr_pres, pres);
313 enerd->term[F_PDISPCORR] = prescorr;
314 enerd->term[F_PRES] += prescorr;
318 /* check whether an 'nst'-style parameter p is a multiple of nst, and
319 set it to be one if not, with a warning. */
320 static void check_nst_param(const gmx::MDLogger &mdlog,
321 const char *desc_nst, int nst,
322 const char *desc_p, int *p)
324 if (*p > 0 && *p % nst != 0)
326 /* Round up to the next multiple of nst */
327 *p = ((*p)/nst + 1)*nst;
328 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
329 "NOTE: %s changes %s to %d", desc_nst, desc_p, *p);
333 void setCurrentLambdasRerun(gmx_int64_t step, const t_lambda *fepvals,
334 const t_trxframe *rerun_fr, const double *lam0,
335 t_state *globalState)
337 GMX_RELEASE_ASSERT(globalState != nullptr, "setCurrentLambdasGlobalRerun should be called with a valid state object");
339 if (rerun_fr->bLambda)
341 if (fepvals->delta_lambda == 0)
343 globalState->lambda[efptFEP] = rerun_fr->lambda;
347 /* find out between which two value of lambda we should be */
348 real frac = step*fepvals->delta_lambda;
349 int fep_state = static_cast<int>(std::floor(frac*fepvals->n_lambda));
350 /* interpolate between this state and the next */
351 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
352 frac = frac*fepvals->n_lambda - fep_state;
353 for (int i = 0; i < efptNR; i++)
355 globalState->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
356 frac*(fepvals->all_lambda[i][fep_state+1] - fepvals->all_lambda[i][fep_state]);
360 else if (rerun_fr->bFepState)
362 globalState->fep_state = rerun_fr->fep_state;
363 for (int i = 0; i < efptNR; i++)
365 globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
370 void setCurrentLambdasLocal(gmx_int64_t step, const t_lambda *fepvals,
371 const double *lam0, t_state *state)
372 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
373 requiring different logic. */
375 if (fepvals->delta_lambda != 0)
377 /* find out between which two value of lambda we should be */
378 real frac = step*fepvals->delta_lambda;
379 if (fepvals->n_lambda > 0)
381 int fep_state = static_cast<int>(std::floor(frac*fepvals->n_lambda));
382 /* interpolate between this state and the next */
383 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
384 frac = frac*fepvals->n_lambda - fep_state;
385 for (int i = 0; i < efptNR; i++)
387 state->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
388 frac*(fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
393 for (int i = 0; i < efptNR; i++)
395 state->lambda[i] = lam0[i] + frac;
401 /* if < 0, fep_state was never defined, and we should not set lambda from the state */
402 if (state->fep_state > -1)
404 for (int i = 0; i < efptNR; i++)
406 state->lambda[i] = fepvals->all_lambda[i][state->fep_state];
412 static void min_zero(int *n, int i)
414 if (i > 0 && (*n == 0 || i < *n))
420 static int lcd4(int i1, int i2, int i3, int i4)
431 gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
434 while (nst > 1 && ((i1 > 0 && i1 % nst != 0) ||
435 (i2 > 0 && i2 % nst != 0) ||
436 (i3 > 0 && i3 % nst != 0) ||
437 (i4 > 0 && i4 % nst != 0)))
445 int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputrec *ir)
447 if (!EI_DYNAMICS(ir->eI))
452 if (nstglobalcomm == -1)
454 // Set up the default behaviour
455 if (!(ir->nstcalcenergy > 0 ||
460 /* The user didn't choose the period for anything
461 important, so we just make sure we can send signals and
462 write output suitably. */
464 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
466 nstglobalcomm = ir->nstenergy;
471 /* The user has made a choice (perhaps implicitly), so we
472 * ensure that we do timely intra-simulation communication
473 * for (possibly) each of the four parts that care.
475 * TODO Does the Verlet scheme (+ DD) need any
476 * communication at nstlist steps? Is the use of nstlist
477 * here a leftover of the twin-range scheme? Can we remove
478 * nstlist when we remove the group scheme?
480 nstglobalcomm = lcd4(ir->nstcalcenergy,
482 ir->etc != etcNO ? ir->nsttcouple : 0,
483 ir->epc != epcNO ? ir->nstpcouple : 0);
488 // Check that the user's choice of mdrun -gcom will work
489 if (ir->nstlist > 0 &&
490 nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
492 nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
493 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
494 "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d",
497 if (ir->nstcalcenergy > 0)
499 check_nst_param(mdlog, "-gcom", nstglobalcomm,
500 "nstcalcenergy", &ir->nstcalcenergy);
502 if (ir->etc != etcNO && ir->nsttcouple > 0)
504 check_nst_param(mdlog, "-gcom", nstglobalcomm,
505 "nsttcouple", &ir->nsttcouple);
507 if (ir->epc != epcNO && ir->nstpcouple > 0)
509 check_nst_param(mdlog, "-gcom", nstglobalcomm,
510 "nstpcouple", &ir->nstpcouple);
513 check_nst_param(mdlog, "-gcom", nstglobalcomm,
514 "nstenergy", &ir->nstenergy);
516 check_nst_param(mdlog, "-gcom", nstglobalcomm,
517 "nstlog", &ir->nstlog);
520 if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
522 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
523 "WARNING: Changing nstcomm from %d to %d",
524 ir->nstcomm, nstglobalcomm);
525 ir->nstcomm = nstglobalcomm;
528 GMX_LOG(mdlog.info).appendTextFormatted(
529 "Intra-simulation communication will occur every %d steps.\n", nstglobalcomm);
530 return nstglobalcomm;
533 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
538 if (MASTER(cr) && *bLastStep)
544 gmx_bcast(sizeof(*fr), fr, cr);
548 *bLastStep = (fr->natoms < 0);
552 // TODO Most of this logic seems to belong in the respective modules
553 void set_state_entries(t_state *state, const t_inputrec *ir)
555 /* The entries in the state in the tpx file might not correspond
556 * with what is needed, so we correct this here.
559 if (ir->efep != efepNO || ir->bExpanded)
561 state->flags |= (1<<estLAMBDA);
562 state->flags |= (1<<estFEPSTATE);
564 state->flags |= (1<<estX);
565 GMX_RELEASE_ASSERT(state->x.size() >= static_cast<unsigned int>(state->natoms), "We should start a run with an initialized state->x");
566 if (EI_DYNAMICS(ir->eI))
568 state->flags |= (1<<estV);
572 if (ir->ePBC != epbcNONE)
574 state->flags |= (1<<estBOX);
575 if (inputrecPreserveShape(ir))
577 state->flags |= (1<<estBOX_REL);
579 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
581 state->flags |= (1<<estBOXV);
582 state->flags |= (1<<estPRES_PREV);
584 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
587 state->flags |= (1<<estNHPRES_XI);
588 state->flags |= (1<<estNHPRES_VXI);
589 state->flags |= (1<<estSVIR_PREV);
590 state->flags |= (1<<estFVIR_PREV);
591 state->flags |= (1<<estVETA);
592 state->flags |= (1<<estVOL0);
594 if (ir->epc == epcBERENDSEN)
596 state->flags |= (1<<estBAROS_INT);
600 if (ir->etc == etcNOSEHOOVER)
602 state->flags |= (1<<estNH_XI);
603 state->flags |= (1<<estNH_VXI);
606 if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
608 state->flags |= (1<<estTHERM_INT);
611 init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
612 init_ekinstate(&state->ekinstate, ir);
616 snew(state->dfhist, 1);
617 init_df_history(state->dfhist, ir->fepvals->n_lambda);