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,2021, 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/utilities.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/mdlib/coupling.h"
54 #include "gromacs/mdlib/gmx_omp_nthreads.h"
55 #include "gromacs/mdlib/simulationsignal.h"
56 #include "gromacs/mdlib/stat.h"
57 #include "gromacs/mdlib/tgroup.h"
58 #include "gromacs/mdlib/update.h"
59 #include "gromacs/mdlib/vcm.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 static void calc_ke_part_normal(gmx::ArrayRef<const gmx::RVec> v,
84 const t_grpopts* opts,
86 gmx_ekindata_t* ekind,
91 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
92 gmx::ArrayRef<t_grp_acc> grpstat = ekind->grpstat;
94 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
95 an option, but not supported now.
96 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
99 /* group velocities are calculated in update_ekindata and
100 * accumulated in acumulate_groups.
101 * Now the partial global and groups ekin.
103 for (g = 0; (g < opts->ngtc); g++)
105 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
108 clear_mat(tcstat[g].ekinf);
109 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
113 clear_mat(tcstat[g].ekinh);
116 ekind->dekindl_old = ekind->dekindl;
117 int nthread = gmx_omp_nthreads_get(emntUpdate);
119 #pragma omp parallel for num_threads(nthread) schedule(static)
120 for (int thread = 0; thread < nthread; thread++)
122 // This OpenMP only loops over arrays and does not call any functions
123 // or memory allocation. It should not be able to throw, so for now
124 // we do not need a try/catch wrapper.
125 int start_t, end_t, n;
133 start_t = ((thread + 0) * md->homenr) / nthread;
134 end_t = ((thread + 1) * md->homenr) / nthread;
136 ekin_sum = ekind->ekin_work[thread];
137 dekindl_sum = ekind->dekindl_work[thread];
139 for (gt = 0; gt < opts->ngtc; gt++)
141 clear_mat(ekin_sum[gt]);
147 for (n = start_t; n < end_t; n++)
157 hm = 0.5 * md->massT[n];
159 for (d = 0; (d < DIM); d++)
161 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
163 for (d = 0; (d < DIM); d++)
165 for (m = 0; (m < DIM); m++)
167 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
168 ekin_sum[gt][m][d] += hm * v_corrt[m] * v_corrt[d];
171 if (md->nMassPerturbed && md->bPerturbed[n])
173 *dekindl_sum += 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
179 for (int thread = 0; thread < nthread; thread++)
181 for (g = 0; g < opts->ngtc; g++)
185 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g], tcstat[g].ekinf);
189 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g], tcstat[g].ekinh);
193 ekind->dekindl += *ekind->dekindl_work[thread];
196 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
199 static void calc_ke_part_visc(const matrix box,
200 gmx::ArrayRef<const gmx::RVec> x,
201 gmx::ArrayRef<const gmx::RVec> v,
202 const t_grpopts* opts,
204 gmx_ekindata_t* ekind,
206 gmx_bool bEkinAveVel)
208 int start = 0, homenr = md->homenr;
209 int g, d, n, m, gt = 0;
212 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
213 t_cos_acc* cosacc = &(ekind->cosacc);
218 for (g = 0; g < opts->ngtc; g++)
220 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
221 clear_mat(ekind->tcstat[g].ekinh);
223 ekind->dekindl_old = ekind->dekindl;
225 fac = 2 * M_PI / box[ZZ][ZZ];
228 for (n = start; n < start + homenr; n++)
234 hm = 0.5 * md->massT[n];
236 /* Note that the times of x and v differ by half a step */
237 /* MRS -- would have to be changed for VV */
238 cosz = std::cos(fac * x[n][ZZ]);
239 /* Calculate the amplitude of the new velocity profile */
240 mvcos += 2 * cosz * md->massT[n] * v[n][XX];
242 copy_rvec(v[n], v_corrt);
243 /* Subtract the profile for the kinetic energy */
244 v_corrt[XX] -= cosz * cosacc->vcos;
245 for (d = 0; (d < DIM); d++)
247 for (m = 0; (m < DIM); m++)
249 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
252 tcstat[gt].ekinf[m][d] += hm * v_corrt[m] * v_corrt[d];
256 tcstat[gt].ekinh[m][d] += hm * v_corrt[m] * v_corrt[d];
260 if (md->nPerturbed && md->bPerturbed[n])
262 /* The minus sign here might be confusing.
263 * The kinetic contribution from dH/dl doesn't come from
264 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
265 * where p are the momenta. The difference is only a minus sign.
267 dekindl -= 0.5 * (md->massB[n] - md->massA[n]) * iprod(v_corrt, v_corrt);
270 ekind->dekindl = dekindl;
271 cosacc->mvcos = mvcos;
273 inc_nrnb(nrnb, eNR_EKIN, homenr);
276 static void calc_ke_part(gmx::ArrayRef<const gmx::RVec> x,
277 gmx::ArrayRef<const gmx::RVec> v,
279 const t_grpopts* opts,
281 gmx_ekindata_t* ekind,
283 gmx_bool bEkinAveVel)
285 if (ekind->cosacc.cos_accel == 0)
287 calc_ke_part_normal(v, opts, md, ekind, nrnb, bEkinAveVel);
291 calc_ke_part_visc(box, x, v, opts, md, ekind, nrnb, bEkinAveVel);
295 /* TODO Specialize this routine into init-time and loop-time versions?
296 e.g. bReadEkin is only true when restoring from checkpoint */
297 void compute_globals(gmx_global_stat* gstat,
299 const t_inputrec* ir,
301 gmx_ekindata_t* ekind,
302 gmx::ArrayRef<const gmx::RVec> x,
303 gmx::ArrayRef<const gmx::RVec> v,
305 const t_mdatoms* mdatoms,
308 gmx_wallcycle_t wcycle,
309 gmx_enerdata_t* enerd,
314 gmx::Constraints* constr,
315 gmx::SimulationSignaller* signalCoordinator,
316 const matrix lastbox,
317 int* totalNumberOfBondedInteractions,
318 gmx_bool* bSumEkinhOld,
321 gmx_bool bEner, bPres, bTemp;
322 gmx_bool bStopCM, bGStat, bReadEkin, bEkinAveVel, bScaleEkin, bConstrain;
323 gmx_bool bCheckNumberOfBondedInteractions;
326 /* translate CGLO flags to gmx_booleans */
327 bStopCM = ((flags & CGLO_STOPCM) != 0);
328 bGStat = ((flags & CGLO_GSTAT) != 0);
329 bReadEkin = ((flags & CGLO_READEKIN) != 0);
330 bScaleEkin = ((flags & CGLO_SCALEEKIN) != 0);
331 bEner = ((flags & CGLO_ENERGY) != 0);
332 bTemp = ((flags & CGLO_TEMPERATURE) != 0);
333 bPres = ((flags & CGLO_PRESSURE) != 0);
334 bConstrain = ((flags & CGLO_CONSTRAINT) != 0);
335 bCheckNumberOfBondedInteractions = ((flags & CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS) != 0);
337 /* we calculate a full state kinetic energy either with full-step velocity verlet
338 or half step where we need the pressure */
340 bEkinAveVel = (ir->eI == eiVV || (ir->eI == eiVVAK && bPres) || bReadEkin);
342 /* in initalization, it sums the shake virial in vv, and to
343 sums ekinh_old in leapfrog (or if we are calculating ekinh_old) for other reasons */
345 /* ########## Kinetic energy ############## */
349 /* Non-equilibrium MD: this is parallellized, but only does communication
350 * when there really is NEMD.
353 if (PAR(cr) && (ekind->bNEMD))
355 accumulate_u(cr, &(ir->opts), ekind);
359 calc_ke_part(x, v, box, &(ir->opts), mdatoms, ekind, nrnb, bEkinAveVel);
363 /* Calculate center of mass velocity if necessary, also parallellized */
366 calc_vcm_grp(*mdatoms, x, v, vcm);
369 if (bTemp || bStopCM || bPres || bEner || bConstrain || bCheckNumberOfBondedInteractions)
373 /* We will not sum ekinh_old,
374 * so signal that we still have to do it.
376 *bSumEkinhOld = TRUE;
380 gmx::ArrayRef<real> signalBuffer = signalCoordinator->getCommunicationBuffer();
383 wallcycle_start(wcycle, ewcMoveE);
392 bStopCM ? vcm : nullptr,
395 totalNumberOfBondedInteractions,
398 wallcycle_stop(wcycle, ewcMoveE);
400 signalCoordinator->finalizeSignals();
401 *bSumEkinhOld = FALSE;
407 /* Calculate the amplitude of the cosine velocity profile */
408 ekind->cosacc.vcos = ekind->cosacc.mvcos / mdatoms->tmass;
413 /* Sum the kinetic energies of the groups & calc temp */
414 /* compute full step kinetic energies if vv, or if vv-avek and we are computing the pressure with inputrecNptTrotter */
415 /* three maincase: VV with AveVel (md-vv), vv with AveEkin (md-vv-avek), leap with AveEkin (md).
416 Leap with AveVel is not supported; it's not clear that it will actually work.
417 bEkinAveVel: If TRUE, we simply multiply ekin by ekinscale to get a full step kinetic energy.
418 If FALSE, we average ekinh_old and ekinh*ekinscale_nhc to get an averaged half step kinetic energy.
420 enerd->term[F_TEMP] = sum_ekin(&(ir->opts), ekind, &dvdl_ekin, bEkinAveVel, bScaleEkin);
421 enerd->dvdl_lin[efptMASS] = static_cast<double>(dvdl_ekin);
423 enerd->term[F_EKIN] = trace(ekind->ekin);
426 /* ########## Now pressure ############## */
427 // TODO: For the VV integrator bConstrain is needed in the conditional. This is confusing, so get rid of this.
428 if (bPres || bConstrain)
430 m_add(force_vir, shake_vir, total_vir);
432 /* Calculate pressure and apply LR correction if PPPM is used.
433 * Use the box from last timestep since we already called update().
436 enerd->term[F_PRES] = calc_pres(fr->pbcType, ir->nwall, lastbox, ekind->ekin, total_vir, pres);
440 static void min_zero(int* n, int i)
442 if (i > 0 && (*n == 0 || i < *n))
448 static int lcd4(int i1, int i2, int i3, int i4)
459 gmx_incons("All 4 inputs for determining nstglobalcomm are <= 0");
463 && ((i1 > 0 && i1 % nst != 0) || (i2 > 0 && i2 % nst != 0) || (i3 > 0 && i3 % nst != 0)
464 || (i4 > 0 && i4 % nst != 0)))
472 int computeGlobalCommunicationPeriod(const t_inputrec* ir)
474 int nstglobalcomm = 10;
476 // Set up the default behaviour
477 if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != etcNO || ir->epc != epcNO))
479 /* The user didn't choose the period for anything
480 important, so we just make sure we can send signals and
481 write output suitably. */
482 if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
484 nstglobalcomm = ir->nstenergy;
489 /* The user has made a choice (perhaps implicitly), so we
490 * ensure that we do timely intra-simulation communication
491 * for (possibly) each of the four parts that care.
493 * TODO Does the Verlet scheme (+ DD) need any
494 * communication at nstlist steps? Is the use of nstlist
495 * here a leftover of the twin-range scheme? Can we remove
496 * nstlist when we remove the group scheme?
498 nstglobalcomm = lcd4(ir->nstcalcenergy,
500 ir->etc != etcNO ? ir->nsttcouple : 0,
501 ir->epc != epcNO ? ir->nstpcouple : 0);
504 return nstglobalcomm;
507 int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, const t_inputrec* ir, const t_commrec* cr)
509 const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
514 .appendTextFormatted("Intra-simulation communication will occur every %d steps.\n",
517 return nstglobalcomm;
520 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
524 if (MASTER(cr) && *bLastStep)
530 gmx_bcast(sizeof(*fr), fr, cr->mpi_comm_mygroup);
534 *bLastStep = (fr->natoms < 0);
537 // TODO Most of this logic seems to belong in the respective modules
538 void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
540 /* The entries in the state in the tpx file might not correspond
541 * with what is needed, so we correct this here.
544 if (ir->efep != efepNO || ir->bExpanded)
546 state->flags |= (1 << estLAMBDA);
547 state->flags |= (1 << estFEPSTATE);
549 state->flags |= (1 << estX);
550 GMX_RELEASE_ASSERT(state->x.size() == state->natoms,
551 "We should start a run with an initialized state->x");
552 if (EI_DYNAMICS(ir->eI))
554 state->flags |= (1 << estV);
558 if (ir->pbcType != PbcType::No)
560 state->flags |= (1 << estBOX);
561 if (inputrecPreserveShape(ir))
563 state->flags |= (1 << estBOX_REL);
565 if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
567 state->flags |= (1 << estBOXV);
568 if (!useModularSimulator)
570 state->flags |= (1 << estPRES_PREV);
573 if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
576 state->flags |= (1 << estNHPRES_XI);
577 state->flags |= (1 << estNHPRES_VXI);
578 state->flags |= (1 << estSVIR_PREV);
579 state->flags |= (1 << estFVIR_PREV);
580 state->flags |= (1 << estVETA);
581 state->flags |= (1 << estVOL0);
583 if (ir->epc == epcBERENDSEN || ir->epc == epcCRESCALE)
585 state->flags |= (1 << estBAROS_INT);
589 if (ir->etc == etcNOSEHOOVER)
591 state->flags |= (1 << estNH_XI);
592 state->flags |= (1 << estNH_VXI);
595 if (ir->etc == etcVRESCALE || ir->etc == etcBERENDSEN)
597 state->flags |= (1 << estTHERM_INT);
600 init_gtc_state(state, state->ngtc, state->nnhpres, ir->opts.nhchainlength); /* allocate the space for nose-hoover chains */
601 init_ekinstate(&state->ekinstate, ir);
605 snew(state->dfhist, 1);
606 init_df_history(state->dfhist, ir->fepvals->n_lambda);
609 if (ir->pull && ir->pull->bSetPbcRefToPrevStepCOM)
611 state->flags |= (1 << estPULLCOMPREVSTEP);