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 The GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "md_support.h"
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/coupling.h"
53 #include "gromacs/mdlib/gmx_omp_nthreads.h"
54 #include "gromacs/mdlib/simulationsignal.h"
55 #include "gromacs/mdlib/stat.h"
56 #include "gromacs/mdlib/tgroup.h"
57 #include "gromacs/mdlib/update.h"
58 #include "gromacs/mdlib/vcm.h"
59 #include "gromacs/mdrunutility/multisim.h"
60 #include "gromacs/mdtypes/commrec.h"
61 #include "gromacs/mdtypes/df_history.h"
62 #include "gromacs/mdtypes/enerdata.h"
63 #include "gromacs/mdtypes/energyhistory.h"
64 #include "gromacs/mdtypes/forcerec.h"
65 #include "gromacs/mdtypes/group.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/mdatom.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/pbc.h"
71 #include "gromacs/pulling/pull.h"
72 #include "gromacs/timing/wallcycle.h"
73 #include "gromacs/topology/mtop_util.h"
74 #include "gromacs/trajectory/trajectoryframe.h"
75 #include "gromacs/utility/arrayref.h"
76 #include "gromacs/utility/cstringutil.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/gmxassert.h"
79 #include "gromacs/utility/logger.h"
80 #include "gromacs/utility/smalloc.h"
81 #include "gromacs/utility/snprintf.h"
83 // TODO move this to multi-sim module
84 bool multisim_int_all_are_equal(const gmx_multisim_t* ms, int64_t value)
86 bool allValuesAreEqual = true;
89 GMX_RELEASE_ASSERT(ms, "Invalid use of multi-simulation pointer");
92 /* send our value to all other master ranks, receive all of theirs */
94 gmx_sumli_sim(ms->nsim, buf, ms);
96 for (int s = 0; s < ms->nsim; s++)
100 allValuesAreEqual = false;
107 return allValuesAreEqual;
110 int multisim_min(const gmx_multisim_t* ms, int nmin, int n)
113 gmx_bool bPos, bEqual;
118 gmx_sumi_sim(ms->nsim, buf, ms);
121 for (s = 0; s < ms->nsim; s++)
123 bPos = bPos && (buf[s] > 0);
124 bEqual = bEqual && (buf[s] == buf[0]);
130 nmin = std::min(nmin, buf[0]);
134 /* Find the least common multiple */
135 for (d = 2; d < nmin; d++)
138 while (s < ms->nsim && d % buf[s] == 0)
144 /* We found the LCM and it is less than nmin */
156 static void calc_ke_part_normal(gmx::ArrayRef<const gmx::RVec> v,
157 const t_grpopts* opts,
159 gmx_ekindata_t* ekind,
161 gmx_bool bEkinAveVel)
164 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
165 gmx::ArrayRef<t_grp_acc> grpstat = ekind->grpstat;
167 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
168 an option, but not supported now.
169 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
172 /* group velocities are calculated in update_ekindata and
173 * accumulated in acumulate_groups.
174 * Now the partial global and groups ekin.
176 for (g = 0; (g < opts->ngtc); g++)
178 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
181 clear_mat(tcstat[g].ekinf);
182 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
186 clear_mat(tcstat[g].ekinh);
189 ekind->dekindl_old = ekind->dekindl;
190 int nthread = gmx_omp_nthreads_get(emntUpdate);
192 #pragma omp parallel for num_threads(nthread) schedule(static)
193 for (int thread = 0; thread < nthread; thread++)
195 // This OpenMP only loops over arrays and does not call any functions
196 // or memory allocation. It should not be able to throw, so for now
197 // we do not need a try/catch wrapper.
198 int start_t, end_t, n;
206 start_t = ((thread + 0) * md->homenr) / nthread;
207 end_t = ((thread + 1) * md->homenr) / nthread;
209 ekin_sum = ekind->ekin_work[thread];
210 dekindl_sum = ekind->dekindl_work[thread];
212 for (gt = 0; gt < opts->ngtc; gt++)
214 clear_mat(ekin_sum[gt]);
220 for (n = start_t; n < end_t; n++)
230 hm = 0.5 * md->massT[n];
232 for (d = 0; (d < DIM); d++)
234 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
236 for (d = 0; (d < DIM); d++)
238 for (m = 0; (m < DIM); m++)
240 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
241 ekin_sum[gt][m][d] += hm * v_corrt[m] * v_corrt[d];
244 if (md->nMassPerturbed && md->bPerturbed[n])
246 *dekindl_sum += 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
252 for (int thread = 0; thread < nthread; thread++)
254 for (g = 0; g < opts->ngtc; g++)
258 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g], tcstat[g].ekinf);
262 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g], tcstat[g].ekinh);
266 ekind->dekindl += *ekind->dekindl_work[thread];
269 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
272 static void calc_ke_part_visc(const matrix box,
273 gmx::ArrayRef<const gmx::RVec> x,
274 gmx::ArrayRef<const gmx::RVec> v,
275 const t_grpopts* opts,
277 gmx_ekindata_t* ekind,
279 gmx_bool bEkinAveVel)
281 int start = 0, homenr = md->homenr;
282 int g, d, n, m, gt = 0;
285 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
286 t_cos_acc* cosacc = &(ekind->cosacc);
291 for (g = 0; g < opts->ngtc; g++)
293 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
294 clear_mat(ekind->tcstat[g].ekinh);
296 ekind->dekindl_old = ekind->dekindl;
298 fac = 2 * M_PI / box[ZZ][ZZ];
301 for (n = start; n < start + homenr; n++)
307 hm = 0.5 * md->massT[n];
309 /* Note that the times of x and v differ by half a step */
310 /* MRS -- would have to be changed for VV */
311 cosz = std::cos(fac * x[n][ZZ]);
312 /* Calculate the amplitude of the new velocity profile */
313 mvcos += 2 * cosz * md->massT[n] * v[n][XX];
315 copy_rvec(v[n], v_corrt);
316 /* Subtract the profile for the kinetic energy */
317 v_corrt[XX] -= cosz * cosacc->vcos;
318 for (d = 0; (d < DIM); d++)
320 for (m = 0; (m < DIM); m++)
322 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
325 tcstat[gt].ekinf[m][d] += hm * v_corrt[m] * v_corrt[d];
329 tcstat[gt].ekinh[m][d] += hm * v_corrt[m] * v_corrt[d];
333 if (md->nPerturbed && md->bPerturbed[n])
335 /* The minus sign here might be confusing.
336 * The kinetic contribution from dH/dl doesn't come from
337 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
338 * where p are the momenta. The difference is only a minus sign.
340 dekindl -= 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
343 ekind->dekindl = dekindl;
344 cosacc->mvcos = mvcos;
346 inc_nrnb(nrnb, eNR_EKIN, homenr);
349 static void calc_ke_part(gmx::ArrayRef<const gmx::RVec> x,
350 gmx::ArrayRef<const gmx::RVec> v,
352 const t_grpopts* opts,
354 gmx_ekindata_t* ekind,
356 gmx_bool bEkinAveVel)
358 if (ekind->cosacc.cos_accel == 0)
360 calc_ke_part_normal(v, opts, md, ekind, nrnb, bEkinAveVel);
364 calc_ke_part_visc(box, x, v, opts, md, ekind, nrnb, bEkinAveVel);
368 /* TODO Specialize this routine into init-time and loop-time versions?
369 e.g. bReadEkin is only true when restoring from checkpoint */
370 void compute_globals(gmx_global_stat* gstat,
372 const t_inputrec* ir,
374 gmx_ekindata_t* ekind,
375 gmx::ArrayRef<const gmx::RVec> x,
376 gmx::ArrayRef<const gmx::RVec> v,
378 const t_mdatoms* mdatoms,
381 gmx_wallcycle_t wcycle,
382 gmx_enerdata_t* enerd,
387 gmx::Constraints* constr,
388 gmx::SimulationSignaller* signalCoordinator,
389 const matrix lastbox,
390 int* totalNumberOfBondedInteractions,
391 gmx_bool* bSumEkinhOld,
394 gmx_bool bEner, bPres, bTemp;
395 gmx_bool bStopCM, bGStat, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
396 gmx_bool bCheckNumberOfBondedInteractions;
399 /* translate CGLO flags to gmx_booleans */
400 bStopCM = ((flags & CGLO_STOPCM) != 0);
401 bGStat = ((flags & CGLO_GSTAT) != 0);
402 bReadEkin = ((flags & CGLO_READEKIN) != 0);
403 bScaleEkin = ((flags & CGLO_SCALEEKIN) != 0);
404 bEner = ((flags & CGLO_ENERGY) != 0);
405 bTemp = ((flags & CGLO_TEMPERATURE) != 0);
406 bPres = ((flags & CGLO_PRESSURE) != 0);
407 bConstrain = ((flags & CGLO_CONSTRAINT) != 0);
408 bCheckNumberOfBondedInteractions = ((flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0);
410 /* we calculate a full state kinetic energy either with full-step velocity verlet
411 or half step where we need the pressure */
413 bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
415 /* in initalization, it sums the shake virial in vv, and to
416 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
418 /* ########## Kinetic energy ############## */
422 /* Non-equilibrium MD: this is parallellized, but only does communication
423 * when there really is NEMD.
426 if (PAR(cr) && (ekind->bNEMD))
428 accumulate_u(cr, &(ir->opts), ekind);
432 calc_ke_part(x, v, box, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
436 /* Calculate center of mass velocity if necessary, also parallellized */
439 calc_vcm_grp(*mdatoms, x, v, vcm);
442 if (bTemp || bStopCM || bPres || bEner || bConstrain || bCheckNumberOfBondedInteractions)
446 /* We will not sum ekinh_old,
447 * so signal that we still have to do it.
449 *bSumEkinhOld = TRUE;
453 gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
456 wallcycle_start(wcycle, ewcMoveE);
457 global_stat(gstat, cr, enerd, force_vir, shake_vir, ir, ekind, constr,
458 bStopCM ? vcm : nullptr, signalBuffer.size(), signalBuffer.data(),
459 totalNumberOfBondedInteractions, *bSumEkinhOld, flags);
460 wallcycle_stop(wcycle, ewcMoveE);
462 signalCoordinator->finalizeSignals();
463 *bSumEkinhOld = FALSE;
469 /* Calculate the amplitude of the cosine velocity profile */
470 ekind->cosacc.vcos = ekind->cosacc.mvcos / mdatoms->tmass;
475 /* Sum the kinetic energies of the groups & calc temp */
476 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
477 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
478 Leap with AveVel is not supported; it's not clear that it will actually work.
479 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
480 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
482 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin, bEkinAveVel, bScaleEkin);
483 enerd->dvdl_lin[efptMASS] = static_cast<double>(dvdl_ekin);
485 enerd->term[F_EKIN] = trace(ekind->ekin);
487 for (auto& dhdl : enerd->dhdlLambda)
489 dhdl += enerd->dvdl_lin[efptMASS];
493 /* ########## Now pressure ############## */
494 // TODO: For the VV integrator bConstrain is needed in the conditional. This is confusing, so get rid of this.
495 if (bPres || bConstrain)
497 m_add(force_vir, shake_vir, total_vir);
499 /* Calculate pressure and apply LR correction if PPPM is used.
500 * Use the box from last timestep since we already called update().
503 enerd->term[F_PRES] = calc_pres(fr->pbcType, ir->nwall, lastbox, ekind->ekin, total_vir, pres);
507 void setCurrentLambdasRerun(int64_t step,
508 const t_lambda* fepvals,
509 const t_trxframe* rerun_fr,
511 t_state* globalState)
513 GMX_RELEASE_ASSERT(globalState != nullptr,
514 "setCurrentLambdasGlobalRerun should be called with a valid state object");
516 if (rerun_fr->bLambda)
518 if (fepvals->delta_lambda == 0)
520 globalState->lambda[efptFEP] = rerun_fr->lambda;
524 /* find out between which two value of lambda we should be */
525 real frac = step * fepvals->delta_lambda;
526 int fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
527 /* interpolate between this state and the next */
528 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
529 frac = frac * fepvals->n_lambda - fep_state;
530 for (int i = 0; i < efptNR; i++)
532 globalState->lambda[i] =
533 lam0[i] + (fepvals->all_lambda[i][fep_state])
534 + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
538 else if (rerun_fr->bFepState)
540 globalState->fep_state = rerun_fr->fep_state;
541 for (int i = 0; i < efptNR; i++)
543 globalState->lambda[i] = fepvals->all_lambda[i][globalState->fep_state];
548 void setCurrentLambdasLocal(const int64_t step,
549 const t_lambda* fepvals,
551 gmx::ArrayRef<real> lambda,
552 const int currentFEPState)
553 /* find the current lambdas. If rerunning, we either read in a state, or a lambda value,
554 requiring different logic. */
556 if (fepvals->delta_lambda != 0)
558 /* find out between which two value of lambda we should be */
559 real frac = step * fepvals->delta_lambda;
560 if (fepvals->n_lambda > 0)
562 int fep_state = static_cast<int>(std::floor(frac * fepvals->n_lambda));
563 /* interpolate between this state and the next */
564 /* this assumes that the initial lambda corresponds to lambda==0, which is verified in grompp */
565 frac = frac * fepvals->n_lambda - fep_state;
566 for (int i = 0; i < efptNR; i++)
568 lambda[i] = lam0[i] + (fepvals->all_lambda[i][fep_state])
569 + frac * (fepvals->all_lambda[i][fep_state + 1] - fepvals->all_lambda[i][fep_state]);
574 for (int i = 0; i < efptNR; i++)
576 lambda[i] = lam0[i] + frac;
582 /* if < 0, fep_state was never defined, and we should not set lambda from the state */
583 if (currentFEPState > -1)
585 for (int i = 0; i < efptNR; i++)
587 lambda[i] = fepvals->all_lambda[i][currentFEPState];
593 static void min_zero(int* n, int i)
595 if (i > 0 && (*n == 0 || i < *n))
601 static int lcd4(int i1, int i2, int i3, int i4)
612 gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
616 && ((i1 > 0 && i1 % nst != 0) || (i2 > 0 && i2 % nst != 0) || (i3 > 0 && i3 % nst != 0)
617 || (i4 > 0 && i4 % nst != 0)))
625 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr)
629 // Set up the default behaviour
630 if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != etcNO || ir->epc != epcNO))
632 /* The user didn't choose the period for anything
633 important, so we just make sure we can send signals and
634 write output suitably. */
636 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
638 nstglobalcomm = ir->nstenergy;
643 /* The user has made a choice (perhaps implicitly), so we
644 * ensure that we do timely intra-simulation communication
645 * for (possibly) each of the four parts that care.
647 * TODO Does the Verlet scheme (+ DD) need any
648 * communication at nstlist steps? Is the use of nstlist
649 * here a leftover of the twin-range scheme? Can we remove
650 * nstlist when we remove the group scheme?
652 nstglobalcomm = lcd4(ir->nstcalcenergy, ir->nstlist, ir->etc != etcNO ? ir->nsttcouple : 0,
653 ir->epc != epcNO ? ir->nstpcouple : 0);
657 // TODO change this behaviour. Instead grompp should print
658 // a (performance) note and mdrun should not change ir.
659 if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
661 GMX_LOG(mdlog.warning)
663 .appendTextFormatted("WARNING: Changing nstcomm from %d to %d", ir->nstcomm, nstglobalcomm);
664 ir->nstcomm = nstglobalcomm;
670 .appendTextFormatted("Intra-simulation communication will occur every %d steps.\n",
673 return nstglobalcomm;
676 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
680 if (MASTER(cr) && *bLastStep)
686 gmx_bcast(sizeof(*fr), fr, cr->mpi_comm_mygroup);
690 *bLastStep = (fr->natoms < 0);
693 // TODO Most of this logic seems to belong in the respective modules
694 void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
696 /* The entries in the state in the tpx file might not correspond
697 * with what is needed, so we correct this here.
700 if (ir->efep != efepNO || ir->bExpanded)
702 state->flags |= (1 << estLAMBDA);
703 state->flags |= (1 << estFEPSTATE);
705 state->flags |= (1 << estX);
706 GMX_RELEASE_ASSERT(state->x.size() == state->natoms,
707 "We should start a run with an initialized state->x");
708 if (EI_DYNAMICS(ir->eI))
710 state->flags |= (1 << estV);
714 if (ir->pbcType != PbcType::No)
716 state->flags |= (1 << estBOX);
717 if (inputrecPreserveShape(ir))
719 state->flags |= (1 << estBOX_REL);
721 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
723 state->flags |= (1 << estBOXV);
724 if (!useModularSimulator)
726 state->flags |= (1 << estPRES_PREV);
729 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
732 state->flags |= (1 << estNHPRES_XI);
733 state->flags |= (1 << estNHPRES_VXI);
734 state->flags |= (1 << estSVIR_PREV);
735 state->flags |= (1 << estFVIR_PREV);
736 state->flags |= (1 << estVETA);
737 state->flags |= (1 << estVOL0);
739 if (ir->epc == epcBERENDSEN)
741 state->flags |= (1 << estBAROS_INT);
745 if (ir->etc == etcNOSEHOOVER)
747 state->flags |= (1 << estNH_XI);
748 state->flags |= (1 << estNH_VXI);
751 if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
753 state->flags |= (1 << estTHERM_INT);
756 init_gtc_state(state, state->ngtc, state->nnhpres,
757 ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
758 init_ekinstate(&state->ekinstate, ir);
762 snew(state->dfhist, 1);
763 init_df_history(state->dfhist, ir->fepvals->n_lambda);
766 if (ir->pull && ir->pull->bSetPbcRefToPrevStepCOM)
768 state->flags |= (1 << estPULLCOMPREVSTEP);