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.
40 #include "md_support.h"
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,
79 bool allValuesAreEqual = true;
82 GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
85 /* send our value to all other master ranks, receive all of theirs */
87 gmx_sumli_sim(ms->nsim, buf, ms);
89 for (int s = 0; s < ms->nsim; s++)
93 allValuesAreEqual = false;
100 return allValuesAreEqual;
103 int multisim_min(const gmx_multisim_t *ms, int nmin, int n)
106 gmx_bool bPos, bEqual;
111 gmx_sumi_sim(ms->nsim, buf, ms);
114 for (s = 0; s < ms->nsim; s++)
116 bPos = bPos && (buf[s] > 0);
117 bEqual = bEqual && (buf[s] == buf[0]);
123 nmin = std::min(nmin, buf[0]);
127 /* Find the least common multiple */
128 for (d = 2; d < nmin; d++)
131 while (s < ms->nsim && d % buf[s] == 0)
137 /* We found the LCM and it is less than 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, struct gmx_constr *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 ############## */
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);
199 calc_ke_part(state, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
203 /* Calculate center of mass velocity if necessary, also parallellized */
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)
214 /* We will not sum ekinh_old,
215 * so signal that we still have to do it.
217 *bSumEkinhOld = TRUE;
222 gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
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))
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 */
249 check_cm_grp(fplog, vcm, ir, 1);
250 /* At initialization, do not pass x with acceleration-correction mode
251 * to avoid (incorrect) correction of the initial coordinates.
253 rvec *xPtr = nullptr;
254 if (vcm->mode == ecmANGULAR || (vcm->mode == ecmLINEAR_ACCELERATION_CORRECTION && !(flags & CGLO_INITIALIZATION)))
256 xPtr = as_rvec_array(state->x.data());
258 do_stopcm_grp(mdatoms->homenr, mdatoms->cVCM,
259 xPtr, as_rvec_array(state->v.data()), *vcm);
260 inc_nrnb(nrnb, eNR_STOPCM, mdatoms->homenr);
265 /* Calculate the amplitude of the cosine velocity profile */
266 ekind->cosacc.vcos = ekind->cosacc.mvcos/mdatoms->tmass;
271 /* Sum the kinetic energies of the groups & calc temp */
272 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
273 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
274 Leap with AveVel is not supported; it's not clear that it will actually work.
275 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
276 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
278 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin,
279 bEkinAveVel, bScaleEkin);
280 enerd->dvdl_lin[efptMASS] = (double) dvdl_ekin;
282 enerd->term[F_EKIN] = trace(ekind->ekin);
285 /* ########## Long range energy information ###### */
287 if (bEner || bPres || bConstrain)
289 calc_dispcorr(ir, fr, box, state->lambda[efptVDW],
290 corr_pres, corr_vir, &prescorr, &enercorr, &dvdlcorr);
295 enerd->term[F_DISPCORR] = enercorr;
296 enerd->term[F_EPOT] += enercorr;
297 enerd->term[F_DVDL_VDW] += dvdlcorr;
300 /* ########## Now pressure ############## */
301 if (bPres || bConstrain)
304 m_add(force_vir, shake_vir, total_vir);
306 /* Calculate pressure and apply LR correction if PPPM is used.
307 * Use the box from last timestep since we already called update().
310 enerd->term[F_PRES] = calc_pres(fr->ePBC, ir->nwall, box, ekind->ekin, total_vir, pres);
312 /* Calculate long range corrections to pressure and energy */
313 /* this adds to enerd->term[F_PRES] and enerd->term[F_ETOT],
314 and computes enerd->term[F_DISPCORR]. Also modifies the
315 total_vir and pres tesors */
317 m_add(total_vir, corr_vir, total_vir);
318 m_add(pres, corr_pres, pres);
319 enerd->term[F_PDISPCORR] = prescorr;
320 enerd->term[F_PRES] += prescorr;
324 /* check whether an 'nst'-style parameter p is a multiple of nst, and
325 set it to be one if not, with a warning. */
326 static void check_nst_param(const gmx::MDLogger &mdlog,
327 const char *desc_nst, int nst,
328 const char *desc_p, int *p)
330 if (*p > 0 && *p % nst != 0)
332 /* Round up to the next multiple of nst */
333 *p = ((*p)/nst + 1)*nst;
334 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
335 "NOTE: %s changes %s to %d", desc_nst, desc_p, *p);
339 void setCurrentLambdasRerun(gmx_int64_t step, const t_lambda *fepvals,
340 const t_trxframe *rerun_fr, const double *lam0,
341 t_state *globalState)
343 GMX_RELEASE_ASSERT(globalState != nullptr, "setCurrentLambdasGlobalRerun should be called with a valid state object");
345 if (rerun_fr->bLambda)
347 if (fepvals->delta_lambda == 0)
349 globalState->lambda[efptFEP] = rerun_fr->lambda;
353 /* find out between which two value of lambda we should be */
354 real frac = step*fepvals->delta_lambda;
355 int fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
356 /* interpolate between this state and the next */
357 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
358 frac = frac*fepvals->n_lambda - fep_state;
359 for (int i = 0; i < efptNR; i++)
361 globalState->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
362 frac*(fepvals->all_lambda[i][fep_state+1] - fepvals->all_lambda[i][fep_state]);
366 else if (rerun_fr->bFepState)
368 globalState->fep_state = rerun_fr->fep_state;
369 for (int i = 0; i < efptNR; i++)
371 globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
376 void setCurrentLambdasLocal(gmx_int64_t step, const t_lambda *fepvals,
377 const double *lam0, t_state *state)
378 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
379 requiring different logic. */
381 if (fepvals->delta_lambda != 0)
383 /* find out between which two value of lambda we should be */
384 real frac = step*fepvals->delta_lambda;
385 if (fepvals->n_lambda > 0)
387 int fep_state = static_cast<int>(floor(frac*fepvals->n_lambda));
388 /* interpolate between this state and the next */
389 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
390 frac = frac*fepvals->n_lambda - fep_state;
391 for (int i = 0; i < efptNR; i++)
393 state->lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state]) +
394 frac*(fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
399 for (int i = 0; i < efptNR; i++)
401 state->lambda[i] = lam0[i] + frac;
407 /* if < 0, fep_state was never defined, and we should not set lambda from the state */
408 if (state->fep_state > -1)
410 for (int i = 0; i < efptNR; i++)
412 state->lambda[i] = fepvals->all_lambda[i][state->fep_state];
418 static void min_zero(int *n, int i)
420 if (i > 0 && (*n == 0 || i < *n))
426 static int lcd4(int i1, int i2, int i3, int i4)
437 gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
440 while (nst > 1 && ((i1 > 0 && i1 % nst != 0) ||
441 (i2 > 0 && i2 % nst != 0) ||
442 (i3 > 0 && i3 % nst != 0) ||
443 (i4 > 0 && i4 % nst != 0)))
451 int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputrec *ir)
453 if (!EI_DYNAMICS(ir->eI))
458 if (nstglobalcomm == -1)
460 // Set up the default behaviour
461 if (!(ir->nstcalcenergy > 0 ||
466 /* The user didn't choose the period for anything
467 important, so we just make sure we can send signals and
468 write output suitably. */
470 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
472 nstglobalcomm = ir->nstenergy;
477 /* The user has made a choice (perhaps implicitly), so we
478 * ensure that we do timely intra-simulation communication
479 * for (possibly) each of the four parts that care.
481 * TODO Does the Verlet scheme (+ DD) need any
482 * communication at nstlist steps? Is the use of nstlist
483 * here a leftover of the twin-range scheme? Can we remove
484 * nstlist when we remove the group scheme?
486 nstglobalcomm = lcd4(ir->nstcalcenergy,
488 ir->etc != etcNO ? ir->nsttcouple : 0,
489 ir->epc != epcNO ? ir->nstpcouple : 0);
494 // Check that the user's choice of mdrun -gcom will work
495 if (ir->nstlist > 0 &&
496 nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
498 nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
499 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
500 "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d",
503 if (ir->nstcalcenergy > 0)
505 check_nst_param(mdlog, "-gcom", nstglobalcomm,
506 "nstcalcenergy", &ir->nstcalcenergy);
508 if (ir->etc != etcNO && ir->nsttcouple > 0)
510 check_nst_param(mdlog, "-gcom", nstglobalcomm,
511 "nsttcouple", &ir->nsttcouple);
513 if (ir->epc != epcNO && ir->nstpcouple > 0)
515 check_nst_param(mdlog, "-gcom", nstglobalcomm,
516 "nstpcouple", &ir->nstpcouple);
519 check_nst_param(mdlog, "-gcom", nstglobalcomm,
520 "nstenergy", &ir->nstenergy);
522 check_nst_param(mdlog, "-gcom", nstglobalcomm,
523 "nstlog", &ir->nstlog);
526 if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
528 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
529 "WARNING: Changing nstcomm from %d to %d",
530 ir->nstcomm, nstglobalcomm);
531 ir->nstcomm = nstglobalcomm;
534 GMX_LOG(mdlog.info).appendTextFormatted(
535 "Intra-simulation communication will occur every %d steps.\n", nstglobalcomm);
536 return nstglobalcomm;
539 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
544 if (MASTER(cr) && *bLastStep)
550 gmx_bcast(sizeof(*fr), fr, cr);
554 *bLastStep = (fr->natoms < 0);
558 // TODO Most of this logic seems to belong in the respective modules
559 void set_state_entries(t_state *state, const t_inputrec *ir)
561 /* The entries in the state in the tpx file might not correspond
562 * with what is needed, so we correct this here.
565 if (ir->efep != efepNO || ir->bExpanded)
567 state->flags |= (1<<estLAMBDA);
568 state->flags |= (1<<estFEPSTATE);
570 state->flags |= (1<<estX);
571 GMX_RELEASE_ASSERT(state->x.size() >= static_cast<unsigned int>(state->natoms), "We should start a run with an initialized state->x");
572 if (EI_DYNAMICS(ir->eI))
574 state->flags |= (1<<estV);
578 if (ir->ePBC != epbcNONE)
580 state->flags |= (1<<estBOX);
581 if (inputrecPreserveShape(ir))
583 state->flags |= (1<<estBOX_REL);
585 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
587 state->flags |= (1<<estBOXV);
588 state->flags |= (1<<estPRES_PREV);
590 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
593 state->flags |= (1<<estNHPRES_XI);
594 state->flags |= (1<<estNHPRES_VXI);
595 state->flags |= (1<<estSVIR_PREV);
596 state->flags |= (1<<estFVIR_PREV);
597 state->flags |= (1<<estVETA);
598 state->flags |= (1<<estVOL0);
600 if (ir->epc == epcBERENDSEN)
602 state->flags |= (1<<estBAROS_INT);
606 if (ir->etc == etcNOSEHOOVER)
608 state->flags |= (1<<estNH_XI);
609 state->flags |= (1<<estNH_VXI);
612 if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
614 state->flags |= (1<<estTHERM_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);
622 snew(state->dfhist, 1);
623 init_df_history(state->dfhist, ir->fepvals->n_lambda);