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/dlbtiming.h"
57 #include "gromacs/domdec/domdec.h"
58 #include "gromacs/domdec/domdec_network.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/domdec/mdsetup.h"
61 #include "gromacs/domdec/partition.h"
62 #include "gromacs/essentialdynamics/edsam.h"
63 #include "gromacs/ewald/pme.h"
64 #include "gromacs/ewald/pme_load_balancing.h"
65 #include "gromacs/fileio/trxio.h"
66 #include "gromacs/gmxlib/network.h"
67 #include "gromacs/gmxlib/nrnb.h"
68 #include "gromacs/gpu_utils/gpu_utils.h"
69 #include "gromacs/listed_forces/manage_threading.h"
70 #include "gromacs/math/functions.h"
71 #include "gromacs/math/utilities.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/math/vectypes.h"
74 #include "gromacs/mdlib/checkpointhandler.h"
75 #include "gromacs/mdlib/compute_io.h"
76 #include "gromacs/mdlib/constr.h"
77 #include "gromacs/mdlib/ebin.h"
78 #include "gromacs/mdlib/enerdata_utils.h"
79 #include "gromacs/mdlib/energyoutput.h"
80 #include "gromacs/mdlib/expanded.h"
81 #include "gromacs/mdlib/force.h"
82 #include "gromacs/mdlib/force_flags.h"
83 #include "gromacs/mdlib/forcerec.h"
84 #include "gromacs/mdlib/md_support.h"
85 #include "gromacs/mdlib/mdatoms.h"
86 #include "gromacs/mdlib/mdoutf.h"
87 #include "gromacs/mdlib/membed.h"
88 #include "gromacs/mdlib/resethandler.h"
89 #include "gromacs/mdlib/sighandler.h"
90 #include "gromacs/mdlib/simulationsignal.h"
91 #include "gromacs/mdlib/stat.h"
92 #include "gromacs/mdlib/stophandler.h"
93 #include "gromacs/mdlib/tgroup.h"
94 #include "gromacs/mdlib/trajectory_writing.h"
95 #include "gromacs/mdlib/update.h"
96 #include "gromacs/mdlib/vcm.h"
97 #include "gromacs/mdlib/vsite.h"
98 #include "gromacs/mdrunutility/printtime.h"
99 #include "gromacs/mdtypes/awh_history.h"
100 #include "gromacs/mdtypes/awh_params.h"
101 #include "gromacs/mdtypes/commrec.h"
102 #include "gromacs/mdtypes/df_history.h"
103 #include "gromacs/mdtypes/enerdata.h"
104 #include "gromacs/mdtypes/energyhistory.h"
105 #include "gromacs/mdtypes/fcdata.h"
106 #include "gromacs/mdtypes/forcerec.h"
107 #include "gromacs/mdtypes/group.h"
108 #include "gromacs/mdtypes/inputrec.h"
109 #include "gromacs/mdtypes/interaction_const.h"
110 #include "gromacs/mdtypes/md_enums.h"
111 #include "gromacs/mdtypes/mdatom.h"
112 #include "gromacs/mdtypes/mdrunoptions.h"
113 #include "gromacs/mdtypes/observableshistory.h"
114 #include "gromacs/mdtypes/state.h"
115 #include "gromacs/mimic/communicator.h"
116 #include "gromacs/mimic/utilities.h"
117 #include "gromacs/pbcutil/mshift.h"
118 #include "gromacs/pbcutil/pbc.h"
119 #include "gromacs/pulling/pull.h"
120 #include "gromacs/timing/wallcycle.h"
121 #include "gromacs/timing/walltime_accounting.h"
122 #include "gromacs/topology/atoms.h"
123 #include "gromacs/topology/idef.h"
124 #include "gromacs/topology/mtop_util.h"
125 #include "gromacs/topology/topology.h"
126 #include "gromacs/trajectory/trajectoryframe.h"
127 #include "gromacs/utility/basedefinitions.h"
128 #include "gromacs/utility/cstringutil.h"
129 #include "gromacs/utility/fatalerror.h"
130 #include "gromacs/utility/logger.h"
131 #include "gromacs/utility/real.h"
132 #include "gromacs/utility/smalloc.h"
134 #include "integrator.h"
135 #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 PaddedVector<gmx::RVec> f {};
152 gmx_global_stat_t gstat;
153 t_graph *graph = nullptr;
154 gmx_shellfc_t *shellfc;
158 /* Domain decomposition could incorrectly miss a bonded
159 interaction, but checking for that requires a global
160 communication stage, which does not otherwise happen in DD
161 code. So we do that alongside the first global energy reduction
162 after a new DD is made. These variables handle whether the
163 check happens, and the result it returns. */
164 bool shouldCheckNumberOfBondedInteractions = false;
165 int totalNumberOfBondedInteractions = -1;
167 SimulationSignals signals;
168 // Most global communnication stages don't propagate mdrun
169 // signals, and will use this object to achieve that.
170 SimulationSignaller nullSignaller(nullptr, nullptr, nullptr, false, false);
174 gmx_fatal(FARGS, "Expanded ensemble not supported by MiMiC.");
178 gmx_fatal(FARGS, "Simulated tempering not supported by MiMiC.");
182 gmx_fatal(FARGS, "AWH not supported by MiMiC.");
184 if (replExParams.exchangeInterval > 0)
186 gmx_fatal(FARGS, "Replica exchange not supported by MiMiC.");
188 if (opt2bSet("-ei", nfile, fnm) || observablesHistory->edsamHistory != nullptr)
190 gmx_fatal(FARGS, "Essential dynamics not supported by MiMiC.");
194 gmx_fatal(FARGS, "Interactive MD not supported by MiMiC.");
198 gmx_fatal(FARGS, "Multiple simulations not supported by MiMiC.");
200 if (std::any_of(ir->opts.annealing, ir->opts.annealing + ir->opts.ngtc,
201 [](int i){return i != eannNO; }))
203 gmx_fatal(FARGS, "Simulated annealing not supported by MiMiC.");
206 /* Settings for rerun */
208 ir->nstcalcenergy = 1;
209 int nstglobalcomm = 1;
210 const bool bNS = true;
212 // Communicator to interact with MiMiC
213 MimicCommunicator mimicCommunicator {};
216 mimicCommunicator.init();
217 mimicCommunicator.sendInitData(top_global, state_global->x);
218 ir->nsteps = mimicCommunicator.getStepNumber();
221 ir->nstxout_compressed = 0;
222 SimulationGroups *groups = &top_global->groups;
223 top_global->intermolecularExclusionGroup = genQmmmIndices(*top_global);
225 initialize_lambdas(fplog, *ir, MASTER(cr), &state_global->fep_state, state_global->lambda, lam0);
228 gmx_mdoutf *outf = init_mdoutf(fplog, nfile, fnm, mdrunOptions, cr, outputProvider, ir, top_global, oenv, wcycle);
229 gmx::EnergyOutput energyOutput;
230 energyOutput.prepare(mdoutf_get_fp_ene(outf), top_global, ir, pull_work, mdoutf_get_fp_dhdl(outf), true);
232 /* Kinetic energy data */
233 std::unique_ptr<gmx_ekindata_t> eKinData = std::make_unique<gmx_ekindata_t>();
234 gmx_ekindata_t *ekind = eKinData.get();
235 init_ekindata(fplog, top_global, &(ir->opts), ekind);
236 /* Copy the cos acceleration to the groups struct */
237 ekind->cosacc.cos_accel = ir->cos_accel;
239 gstat = global_stat_init(ir);
241 /* Check for polarizable models and flexible constraints */
242 shellfc = init_shell_flexcon(fplog,
243 top_global, constr ? constr->numFlexibleConstraints() : 0,
244 ir->nstcalcenergy, DOMAINDECOMP(cr));
247 double io = compute_io(ir, top_global->natoms, *groups, energyOutput.numEnergyTerms(), 1);
248 if ((io > 2000) && MASTER(cr))
251 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
256 // Local state only becomes valid now.
257 std::unique_ptr<t_state> stateInstance;
260 if (DOMAINDECOMP(cr))
262 dd_init_local_top(*top_global, &top);
264 stateInstance = std::make_unique<t_state>();
265 state = stateInstance.get();
266 dd_init_local_state(cr->dd, state_global, state);
268 /* Distribute the charge groups over the nodes from the master node */
269 dd_partition_system(fplog, mdlog, ir->init_step, cr, TRUE, 1,
270 state_global, *top_global, ir, imdSession,
272 state, &f, mdAtoms, &top, fr,
274 nrnb, nullptr, FALSE);
275 shouldCheckNumberOfBondedInteractions = true;
276 gmx_bcast(sizeof(ir->nsteps), &ir->nsteps, cr);
280 state_change_natoms(state_global, state_global->natoms);
281 /* We need to allocate one element extra, since we might use
282 * (unaligned) 4-wide SIMD loads to access rvec entries.
284 f.resizeWithPadding(state_global->natoms);
285 /* Copy the pointer to the global state */
286 state = state_global;
288 mdAlgorithmsSetupAtomData(cr, ir, *top_global, &top, fr,
289 &graph, mdAtoms, constr, vsite, shellfc);
292 auto mdatoms = mdAtoms->mdatoms();
294 // NOTE: The global state is no longer used at this point.
295 // But state_global is still used as temporary storage space for writing
296 // the global state to file and potentially for replica exchange.
297 // (Global topology should persist.)
299 update_mdatoms(mdatoms, state->lambda[efptMASS]);
301 if (ir->efep != efepNO && ir->fepvals->nstdhdl != 0)
303 doFreeEnergyPerturbation = true;
307 int cglo_flags = (CGLO_INITIALIZATION | CGLO_GSTAT |
308 (shouldCheckNumberOfBondedInteractions ?
309 CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
310 bool bSumEkinhOld = false;
311 t_vcm *vcm = nullptr;
312 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
313 nullptr, enerd, force_vir, shake_vir, total_vir, pres, mu_tot,
314 constr, &nullSignaller, state->box,
315 &totalNumberOfBondedInteractions, &bSumEkinhOld, cglo_flags);
317 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
318 top_global, &top, state,
319 &shouldCheckNumberOfBondedInteractions);
323 fprintf(stderr, "starting MiMiC MD run '%s'\n\n",
324 *(top_global->name));
325 if (mdrunOptions.verbose)
327 fprintf(stderr, "Calculated time to finish depends on nsteps from "
328 "run input file,\nwhich may not correspond to the time "
329 "needed to process input trajectory.\n\n");
331 fprintf(fplog, "\n");
334 walltime_accounting_start_time(walltime_accounting);
335 wallcycle_start(wcycle, ewcRUN);
336 print_start(fplog, cr, walltime_accounting, "mdrun");
338 /***********************************************************
342 ************************************************************/
346 GMX_LOG(mdlog.info).asParagraph().
347 appendText("Simulations has constraints. Constraints will "
348 "be handled by CPMD.");
351 GMX_LOG(mdlog.info).asParagraph().
352 appendText("MiMiC does not report kinetic energy, total energy, temperature, virial and pressure.");
354 auto stopHandler = stopHandlerBuilder->getStopHandlerMD(
355 compat::not_null<SimulationSignal*>(&signals[eglsSTOPCOND]), false,
356 MASTER(cr), ir->nstlist, mdrunOptions.reproducible, nstglobalcomm,
357 mdrunOptions.maximumHoursToRun, ir->nstlist == 0, fplog, step, bNS, walltime_accounting);
359 // we don't do counter resetting in rerun - finish will always be valid
360 walltime_accounting_set_valid_finish(walltime_accounting);
362 const DDBalanceRegionHandler ddBalanceRegionHandler(cr);
364 step = ir->init_step;
367 /* and stop now if we should */
368 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel > ir->nsteps));
371 isLastStep = (isLastStep || (ir->nsteps >= 0 && step_rel == ir->nsteps));
372 wallcycle_start(wcycle, ewcSTEP);
378 mimicCommunicator.getCoords(&state_global->x, state_global->natoms);
381 if (ir->efep != efepNO)
383 setCurrentLambdasLocal(step, ir->fepvals, lam0, state);
388 const bool constructVsites = ((vsite != nullptr) && mdrunOptions.rerunConstructVsites);
389 if (constructVsites && DOMAINDECOMP(cr))
391 gmx_fatal(FARGS, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
392 "use a single rank");
396 if (DOMAINDECOMP(cr))
398 /* Repartition the domain decomposition */
399 const bool bMasterState = true;
400 dd_partition_system(fplog, mdlog, step, cr,
401 bMasterState, nstglobalcomm,
402 state_global, *top_global, ir, imdSession,
404 state, &f, mdAtoms, &top, fr,
407 mdrunOptions.verbose);
408 shouldCheckNumberOfBondedInteractions = true;
413 print_ebin_header(fplog, step, t); /* can we improve the information printed here? */
416 if (ir->efep != efepNO)
418 update_mdatoms(mdatoms, state->lambda[efptMASS]);
421 force_flags = (GMX_FORCE_STATECHANGED |
422 GMX_FORCE_DYNAMICBOX |
423 GMX_FORCE_ALLFORCES |
424 GMX_FORCE_VIRIAL | // TODO: Get rid of this once #2649 is solved
426 (doFreeEnergyPerturbation ? GMX_FORCE_DHDL : 0));
430 /* Now is the time to relax the shells */
431 relax_shell_flexcon(fplog, cr, ms, mdrunOptions.verbose,
432 enforcedRotation, step,
433 ir, imdSession, pull_work, bNS, force_flags, &top,
435 state, f.arrayRefWithPadding(), force_vir, mdatoms,
437 shellfc, fr, ppForceWorkload, t, mu_tot,
439 ddBalanceRegionHandler);
443 /* The coordinates (x) are shifted (to get whole molecules)
445 * This is parallellized as well, and does communication too.
446 * Check comments in sim_util.c
449 gmx_edsam *ed = nullptr;
450 do_force(fplog, cr, ms, ir, awh, enforcedRotation, imdSession,
452 step, nrnb, wcycle, &top,
453 state->box, state->x.arrayRefWithPadding(), &state->hist,
454 f.arrayRefWithPadding(), force_vir, mdatoms, enerd, fcd,
455 state->lambda, graph,
456 fr, ppForceWorkload, vsite, mu_tot, t, ed,
457 GMX_FORCE_NS | force_flags,
458 ddBalanceRegionHandler);
461 /* Now we have the energies and forces corresponding to the
462 * coordinates at time t.
465 const bool isCheckpointingStep = false;
466 const bool doRerun = false;
467 const bool bSumEkinhOld = false;
468 do_md_trajectory_writing(fplog, cr, nfile, fnm, step, step_rel, t,
469 ir, state, state_global, observablesHistory,
471 outf, energyOutput, ekind, f,
472 isCheckpointingStep, doRerun, isLastStep,
473 mdrunOptions.writeConfout,
477 stopHandler->setSignal();
481 /* Need to unshift here */
482 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
485 if (vsite != nullptr)
487 wallcycle_start(wcycle, ewcVSITECONSTR);
488 if (graph != nullptr)
490 shift_self(graph, state->box, as_rvec_array(state->x.data()));
492 construct_vsites(vsite, as_rvec_array(state->x.data()), ir->delta_t, as_rvec_array(state->v.data()),
493 top.idef.iparams, top.idef.il,
494 fr->ePBC, fr->bMolPBC, cr, state->box);
496 if (graph != nullptr)
498 unshift_self(graph, state->box, as_rvec_array(state->x.data()));
500 wallcycle_stop(wcycle, ewcVSITECONSTR);
504 const bool doInterSimSignal = false;
505 const bool doIntraSimSignal = true;
506 bool bSumEkinhOld = false;
507 t_vcm *vcm = nullptr;
508 SimulationSignaller signaller(&signals, cr, ms, doInterSimSignal, doIntraSimSignal);
510 compute_globals(fplog, gstat, cr, ir, fr, ekind, state, mdatoms, nrnb, vcm,
511 wcycle, enerd, nullptr, nullptr, nullptr, nullptr, mu_tot,
514 &totalNumberOfBondedInteractions, &bSumEkinhOld,
515 CGLO_GSTAT | CGLO_ENERGY
516 | (shouldCheckNumberOfBondedInteractions ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0)
518 checkNumberOfBondedInteractions(mdlog, cr, totalNumberOfBondedInteractions,
519 top_global, &top, state,
520 &shouldCheckNumberOfBondedInteractions);
524 gmx::HostVector<gmx::RVec> fglobal(top_global->natoms);
525 gmx::ArrayRef<gmx::RVec> ftemp;
526 gmx::ArrayRef<const gmx::RVec> flocal = gmx::makeArrayRef(f);
527 if (DOMAINDECOMP(cr))
529 ftemp = gmx::makeArrayRef(fglobal);
530 dd_collect_vec(cr->dd, state, flocal, ftemp);
534 ftemp = gmx::makeArrayRef(f);
539 mimicCommunicator.sendEnergies(enerd->term[F_EPOT]);
540 mimicCommunicator.sendForces(ftemp, state_global->natoms);
546 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
547 the virial that should probably be addressed eventually. state->veta has better properies,
548 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
549 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
551 if (ir->efep != efepNO)
553 /* Sum up the foreign energy and dhdl terms for md and sd.
554 Currently done every step so that dhdl is correct in the .edr */
555 sum_dhdl(enerd, state->lambda, ir->fepvals);
561 const bool bCalcEnerStep = true;
562 energyOutput.addDataAtEnergyStep(doFreeEnergyPerturbation, bCalcEnerStep,
563 t, mdatoms->tmass, enerd, state,
564 ir->fepvals, ir->expandedvals, state->box,
565 shake_vir, force_vir, total_vir, pres,
566 ekind, mu_tot, constr);
568 const bool do_ene = true;
569 const bool do_log = true;
571 const bool do_dr = ir->nstdisreout != 0;
572 const bool do_or = ir->nstorireout != 0;
574 energyOutput.printStepToEnergyFile(mdoutf_get_fp_ene(outf), do_ene, do_dr, do_or,
575 do_log ? fplog : nullptr,
577 eprNORMAL, fcd, groups, &(ir->opts), awh);
579 if (do_per_step(step, ir->nstlog))
581 if (fflush(fplog) != 0)
583 gmx_fatal(FARGS, "Cannot flush logfile - maybe you are out of disk space?");
588 /* Print the remaining wall clock time for the run */
589 if (isMasterSimMasterRank(ms, cr) &&
590 (mdrunOptions.verbose || gmx_got_usr_signal()))
594 fprintf(stderr, "\n");
596 print_time(stderr, walltime_accounting, step, ir, cr);
599 cycles = wallcycle_stop(wcycle, ewcSTEP);
600 if (DOMAINDECOMP(cr) && wcycle)
602 dd_cycles_add(cr->dd, cycles, ddCyclStep);
605 /* increase the MD step number */
609 /* End of main MD loop */
611 /* Closing TNG files can include compressing data. Therefore it is good to do that
612 * before stopping the time measurements. */
613 mdoutf_tng_close(outf);
615 /* Stop measuring walltime */
616 walltime_accounting_end_time(walltime_accounting);
620 mimicCommunicator.finalize();
623 if (!thisRankHasDuty(cr, DUTY_PME))
625 /* Tell the PME only node to finish */
626 gmx_pme_send_finish(cr);
631 done_shellfc(fplog, shellfc, step_rel);
633 walltime_accounting_set_nsteps_done(walltime_accounting, step_rel);