2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * \brief Declares the loop for MiMiC QM/MM
40 * \author Viacheslav Bolnykh <v.bolnykh@hpc-leap.eu>
41 * \ingroup module_mdrun
53 #include "gromacs/awh/awh.h"
54 #include "gromacs/commandline/filenm.h"
55 #include "gromacs/domdec/collect.h"
56 #include "gromacs/domdec/domdec.h"
57 #include "gromacs/domdec/domdec_network.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/partition.h"
60 #include "gromacs/essentialdynamics/edsam.h"
61 #include "gromacs/ewald/pme.h"
62 #include "gromacs/ewald/pme_load_balancing.h"
63 #include "gromacs/fileio/trxio.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/gmxlib/nrnb.h"
66 #include "gromacs/gpu_utils/gpu_utils.h"
67 #include "gromacs/imd/imd.h"
68 #include "gromacs/listed_forces/manage_threading.h"
69 #include "gromacs/math/functions.h"
70 #include "gromacs/math/utilities.h"
71 #include "gromacs/math/vec.h"
72 #include "gromacs/math/vectypes.h"
73 #include "gromacs/mdlib/checkpointhandler.h"
74 #include "gromacs/mdlib/compute_io.h"
75 #include "gromacs/mdlib/constr.h"
76 #include "gromacs/mdlib/ebin.h"
77 #include "gromacs/mdlib/energyoutput.h"
78 #include "gromacs/mdlib/expanded.h"
79 #include "gromacs/mdlib/force.h"
80 #include "gromacs/mdlib/force_flags.h"
81 #include "gromacs/mdlib/forcerec.h"
82 #include "gromacs/mdlib/md_support.h"
83 #include "gromacs/mdlib/mdatoms.h"
84 #include "gromacs/mdlib/mdoutf.h"
85 #include "gromacs/mdlib/mdrun.h"
86 #include "gromacs/mdlib/mdsetup.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/ns.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/shellfc.h"
91 #include "gromacs/mdlib/sighandler.h"
92 #include "gromacs/mdlib/sim_util.h"
93 #include "gromacs/mdlib/simulationsignal.h"
94 #include "gromacs/mdlib/stophandler.h"
95 #include "gromacs/mdlib/tgroup.h"
96 #include "gromacs/mdlib/trajectory_writing.h"
97 #include "gromacs/mdlib/update.h"
98 #include "gromacs/mdlib/vcm.h"
99 #include "gromacs/mdlib/vsite.h"
100 #include "gromacs/mdtypes/awh_history.h"
101 #include "gromacs/mdtypes/awh_params.h"
102 #include "gromacs/mdtypes/commrec.h"
103 #include "gromacs/mdtypes/df_history.h"
104 #include "gromacs/mdtypes/enerdata.h"
105 #include "gromacs/mdtypes/energyhistory.h"
106 #include "gromacs/mdtypes/fcdata.h"
107 #include "gromacs/mdtypes/forcerec.h"
108 #include "gromacs/mdtypes/group.h"
109 #include "gromacs/mdtypes/inputrec.h"
110 #include "gromacs/mdtypes/interaction_const.h"
111 #include "gromacs/mdtypes/md_enums.h"
112 #include "gromacs/mdtypes/mdatom.h"
113 #include "gromacs/mdtypes/observableshistory.h"
114 #include "gromacs/mdtypes/state.h"
115 #include "gromacs/mimic/MimicCommunicator.h"
116 #include "gromacs/mimic/MimicUtils.h"
117 #include "gromacs/pbcutil/mshift.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/swap/swapcoords.h"
121 #include "gromacs/timing/wallcycle.h"
122 #include "gromacs/timing/walltime_accounting.h"
123 #include "gromacs/topology/atoms.h"
124 #include "gromacs/topology/idef.h"
125 #include "gromacs/topology/mtop_util.h"
126 #include "gromacs/topology/topology.h"
127 #include "gromacs/trajectory/trajectoryframe.h"
128 #include "gromacs/utility/basedefinitions.h"
129 #include "gromacs/utility/cstringutil.h"
130 #include "gromacs/utility/fatalerror.h"
131 #include "gromacs/utility/logger.h"
132 #include "gromacs/utility/real.h"
133 #include "gromacs/utility/smalloc.h"
135 #include "integrator.h"
136 #include "replicaexchange.h"
138 using gmx::SimulationSignaller;
140 void gmx::Integrator::do_mimic()
142 t_inputrec *ir = inputrec;
143 int64_t step, step_rel;
144 double t, lam0[efptNR];
145 bool isLastStep = false;
146 bool doFreeEnergyPerturbation = false;
148 tensor force_vir, shake_vir, total_vir, pres;
151 gmx_enerdata_t *enerd;
152 PaddedVector<gmx::RVec> f {};
153 gmx_global_stat_t gstat;
154 t_graph *graph = nullptr;
155 gmx_groups_t *groups;
156 gmx_shellfc_t *shellfc;
160 /* Domain decomposition could incorrectly miss a bonded
161 interaction, but checking for that requires a global
162 communication stage, which does not otherwise happen in DD
163 code. So we do that alongside the first global energy reduction
164 after a new DD is made. These variables handle whether the
165 check happens, and the result it returns. */
166 bool shouldCheckNumberOfBondedInteractions = false;
167 int totalNumberOfBondedInteractions = -1;
169 SimulationSignals signals;
170 // Most global communnication stages don't propagate mdrun
171 // signals, and will use this object to achieve that.
172 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
176 gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
180 gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
184 gmx_fatal(FARGS, "AWH not supported by MiMiC.");
186 if (replExParams.exchangeInterval > 0)
188 gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
190 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
192 gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
196 gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
200 gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
202 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
203 [](int i){return i != eannNO; }))
205 gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
208 /* Settings for rerun */
210 ir->nstcalcenergy = 1;
211 int nstglobalcomm = 1;
212 const bool bNS = true;
214 // Communicator to interact with MiMiC
215 MimicCommunicator mimicCommunicator {};
218 mimicCommunicator.init();
219 mimicCommunicator.sendInitData(top_global, state_global->x);
220 ir->nsteps = mimicCommunicator.getStepNumber();
223 ir->nstxout_compressed = 0;
224 groups = &top_global->groups;
225 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
227 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
230 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle);
231 gmx::EnergyOutput energyOutput;
232 energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, ir, mdoutf_get_fp_dhdl(outf), true);
234 /* Energy terms and groups */
236 init_enerdata(top_global->groups.grps[egcENER].nr, ir->fepvals->n_lambda,
239 /* Kinetic energy data */
240 std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
241 gmx_ekindata_t *ekind = eKinData.get();
242 init_ekindata(fplog, top_global, &(ir->opts), ekind);
243 /* Copy the cos acceleration to the groups struct */
244 ekind->cosacc.cos_accel = ir->cos_accel;
246 gstat = global_stat_init(ir);
248 /* Check for polarizable models and flexible constraints */
249 shellfc = init_shell_flexcon(fplog,
250 top_global, constr ? constr->numFlexibleConstraints() : 0,
251 ir->nstcalcenergy, DOMAINDECOMP(cr));
254 double io = compute_io(ir, top_global->natoms, groups, energyOutput.numEnergyTerms(), 1);
255 if ((io > 2000) && MASTER(cr))
258 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
263 // Local state only becomes valid now.
264 std::unique_ptr<t_state> stateInstance;
267 if (DOMAINDECOMP(cr))
269 dd_init_local_top(*top_global, &top);
271 stateInstance = std::make_unique<t_state>();
272 state = stateInstance.get();
273 dd_init_local_state(cr->dd, state_global, state);
275 /* Distribute the charge groups over the nodes from the master node */
276 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
277 state_global, *top_global, ir,
278 state, &f, mdAtoms, &top, fr,
280 nrnb, nullptr, FALSE);
281 shouldCheckNumberOfBondedInteractions = true;
282 gmx_bcast(sizeof(ir->nsteps), &ir->nsteps, cr);
286 state_change_natoms(state_global, state_global->natoms);
287 /* We need to allocate one element extra, since we might use
288 * (unaligned) 4-wide SIMD loads to access rvec entries.
290 f.resizeWithPadding(state_global->natoms);
291 /* Copy the pointer to the global state */
292 state = state_global;
294 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
295 &graph, mdAtoms, constr, vsite, shellfc);
298 auto mdatoms = mdAtoms->mdatoms();
300 // NOTE: The global state is no longer used at this point.
301 // But state_global is still used as temporary storage space for writing
302 // the global state to file and potentially for replica exchange.
303 // (Global topology should persist.)
305 update_mdatoms(mdatoms, state->lambda[efptMASS]);
307 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
309 doFreeEnergyPerturbation = true;
313 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
314 (shouldCheckNumberOfBondedInteractions ?
315 CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
316 bool bSumEkinhOld = false;
317 t_vcm *vcm = nullptr;
318 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
319 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
320 constr, &nullSignaller, state->box,
321 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
323 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
324 top_global, &top, state,
325 &shouldCheckNumberOfBondedInteractions);
329 fprintf(stderr, "starting MiMiC MD run '%s'\n\n",
330 *(top_global->name));
331 if (mdrunOptions.verbose)
333 fprintf(stderr, "Calculated time to finish depends on nsteps from "
334 "run input file,\nwhich may not correspond to the time "
335 "needed to process input trajectory.\n\n");
337 fprintf(fplog, "\n");
340 walltime_accounting_start_time(walltime_accounting);
341 wallcycle_start(wcycle, ewcRUN);
342 print_start(fplog, cr, walltime_accounting, "mdrun");
344 /***********************************************************
348 ************************************************************/
352 GMX_LOG(mdlog.info).asParagraph().
353 appendText("Simulations has constraints. Constraints will "
354 "be handled by CPMD.");
357 GMX_LOG(mdlog.info).asParagraph().
358 appendText("MiMiC does not report kinetic energy, total energy, temperature, virial and pressure.");
360 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
361 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
362 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
363 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
365 // we don't do counter resetting in rerun - finish will always be valid
366 walltime_accounting_set_valid_finish(walltime_accounting);
368 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion =
370 DdOpenBalanceRegionBeforeForceComputation::yes :
371 DdOpenBalanceRegionBeforeForceComputation::no);
372 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion =
374 DdCloseBalanceRegionAfterForceComputation::yes :
375 DdCloseBalanceRegionAfterForceComputation::no);
377 step = ir->init_step;
380 /* and stop now if we should */
381 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
384 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
385 wallcycle_start(wcycle, ewcSTEP);
391 mimicCommunicator.getCoords(&state_global->x, state_global->natoms);
394 if (ir->efep != efepNO)
396 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
401 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
402 if (constructVsites && DOMAINDECOMP(cr))
404 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
405 "use a single rank");
409 if (DOMAINDECOMP(cr))
411 /* Repartition the domain decomposition */
412 const bool bMasterState = true;
413 dd_partition_system(fplog, mdlog, step, cr,
414 bMasterState, nstglobalcomm,
415 state_global, *top_global, ir,
416 state, &f, mdAtoms, &top, fr,
419 mdrunOptions.verbose);
420 shouldCheckNumberOfBondedInteractions = true;
425 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
428 if (ir->efep != efepNO)
430 update_mdatoms(mdatoms, state->lambda[efptMASS]);
433 force_flags = (GMX_FORCE_STATECHANGED |
434 GMX_FORCE_DYNAMICBOX |
435 GMX_FORCE_ALLFORCES |
436 GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
438 (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
442 /* Now is the time to relax the shells */
443 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
444 enforcedRotation, step,
445 ir, bNS, force_flags, &top,
447 state, f.arrayRefWithPadding(), force_vir, mdatoms,
448 nrnb, wcycle, graph, groups,
449 shellfc, fr, ppForceWorkload, t, mu_tot,
451 ddOpenBalanceRegion, ddCloseBalanceRegion);
455 /* The coordinates (x) are shifted (to get whole molecules)
457 * This is parallellized as well, and does communication too.
458 * Check comments in sim_util.c
461 gmx_edsam *ed = nullptr;
462 do_force(fplog, cr, ms, ir, awh, enforcedRotation,
463 step, nrnb, wcycle, &top, groups,
464 state->box, state->x.arrayRefWithPadding(), &state->hist,
465 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
466 state->lambda, graph,
467 fr, ppForceWorkload, vsite, mu_tot, t, ed,
468 GMX_FORCE_NS | force_flags,
469 ddOpenBalanceRegion, ddCloseBalanceRegion);
472 /* Now we have the energies and forces corresponding to the
473 * coordinates at time t.
476 const bool isCheckpointingStep = false;
477 const bool doRerun = false;
478 const bool bSumEkinhOld = false;
479 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
480 ir, state, state_global, observablesHistory,
482 outf, energyOutput, ekind, f,
483 isCheckpointingStep, doRerun, isLastStep,
484 mdrunOptions.writeConfout,
488 stopHandler->setSignal();
492 /* Need to unshift here */
493 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
496 if (vsite != nullptr)
498 wallcycle_start(wcycle, ewcVSITECONSTR);
499 if (graph != nullptr)
501 shift_self(graph, state->box, as_rvec_array(state->x.data()));
503 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
504 top.idef.iparams, top.idef.il,
505 fr->ePBC, fr->bMolPBC, cr, state->box);
507 if (graph != nullptr)
509 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
511 wallcycle_stop(wcycle, ewcVSITECONSTR);
515 const bool doInterSimSignal = false;
516 const bool doIntraSimSignal = true;
517 bool bSumEkinhOld = false;
518 t_vcm *vcm = nullptr;
519 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
521 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
522 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
525 &totalNumberOfBondedInteractions, &bSumEkinhOld,
526 CGLO_GSTAT | CGLO_ENERGY
527 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
529 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
530 top_global, &top, state,
531 &shouldCheckNumberOfBondedInteractions);
535 gmx::HostVector<gmx::RVec> fglobal(top_global->natoms);
536 gmx::ArrayRef<gmx::RVec> ftemp;
537 gmx::ArrayRef<const gmx::RVec> flocal = gmx::makeArrayRef(f);
538 if (DOMAINDECOMP(cr))
540 ftemp = gmx::makeArrayRef(fglobal);
541 dd_collect_vec(cr->dd, state, flocal, ftemp);
545 ftemp = gmx::makeArrayRef(f);
550 mimicCommunicator.sendEnergies(enerd->term[F_EPOT]);
551 mimicCommunicator.sendForces(ftemp, state_global->natoms);
557 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
558 the virial that should probably be addressed eventually. state->veta has better properies,
559 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
560 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
562 if (ir->efep != efepNO)
564 /* Sum up the foreign energy and dhdl terms for md and sd.
565 Currently done every step so that dhdl is correct in the .edr */
566 sum_dhdl(enerd, state->lambda, ir->fepvals);
572 const bool bCalcEnerStep = true;
573 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
574 t, mdatoms->tmass, enerd, state,
575 ir->fepvals, ir->expandedvals, state->box,
576 shake_vir, force_vir, total_vir, pres,
577 ekind, mu_tot, constr);
579 const bool do_ene = true;
580 const bool do_log = true;
582 const bool do_dr = ir->nstdisreout != 0;
583 const bool do_or = ir->nstorireout != 0;
585 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
586 do_log ? fplog : nullptr,
588 eprNORMAL, fcd, groups, &(ir->opts), awh);
590 if (do_per_step(step, ir->nstlog))
592 if (fflush(fplog) != 0)
594 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
599 /* Print the remaining wall clock time for the run */
600 if (isMasterSimMasterRank(ms, cr) &&
601 (mdrunOptions.verbose || gmx_got_usr_signal()))
605 fprintf(stderr, "\n");
607 print_time(stderr, walltime_accounting, step, ir, cr);
610 cycles = wallcycle_stop(wcycle, ewcSTEP);
611 if (DOMAINDECOMP(cr) && wcycle)
613 dd_cycles_add(cr->dd, cycles, ddCyclStep);
616 /* increase the MD step number */
620 /* End of main MD loop */
622 /* Closing TNG files can include compressing data. Therefore it is good to do that
623 * before stopping the time measurements. */
624 mdoutf_tng_close(outf);
626 /* Stop measuring walltime */
627 walltime_accounting_end_time(walltime_accounting);
631 mimicCommunicator.finalize();
634 if (!thisRankHasDuty(cr, DUTY_PME))
636 /* Tell the PME only node to finish */
637 gmx_pme_send_finish(cr);
642 done_shellfc(fplog, shellfc, step_rel);
644 // Clean up swapcoords
645 if (ir->eSwapCoords != eswapNO)
647 finish_swapcoords(ir->swap);
650 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);
652 destroy_enerdata(enerd);