2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020, 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_load_balancing.h"
64 #include "gromacs/ewald/pme_pp.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/listed_forces.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/handlerestart.h"
99 #include "gromacs/mdrunutility/multisim.h"
100 #include "gromacs/mdrunutility/printtime.h"
101 #include "gromacs/mdtypes/awh_history.h"
102 #include "gromacs/mdtypes/awh_params.h"
103 #include "gromacs/mdtypes/commrec.h"
104 #include "gromacs/mdtypes/df_history.h"
105 #include "gromacs/mdtypes/enerdata.h"
106 #include "gromacs/mdtypes/energyhistory.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/mdrunoptions.h"
114 #include "gromacs/mdtypes/observableshistory.h"
115 #include "gromacs/mdtypes/state.h"
116 #include "gromacs/mimic/communicator.h"
117 #include "gromacs/mimic/utilities.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 "legacysimulator.h"
135 #include "replicaexchange.h"
138 using gmx::SimulationSignaller
;
140 void gmx::LegacySimulator::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;
147 unsigned int force_flags
;
148 tensor force_vir
, shake_vir
, total_vir
, pres
;
150 PaddedHostVector
<gmx::RVec
> f
{};
151 gmx_global_stat_t gstat
;
152 gmx_shellfc_t
* shellfc
;
156 /* Domain decomposition could incorrectly miss a bonded
157 interaction, but checking for that requires a global
158 communication stage, which does not otherwise happen in DD
159 code. So we do that alongside the first global energy reduction
160 after a new DD is made. These variables handle whether the
161 check happens, and the result it returns. */
162 bool shouldCheckNumberOfBondedInteractions
= false;
163 int totalNumberOfBondedInteractions
= -1;
165 SimulationSignals signals
;
166 // Most global communnication stages don't propagate mdrun
167 // signals, and will use this object to achieve that.
168 SimulationSignaller
nullSignaller(nullptr, nullptr, nullptr, false, false);
172 gmx_fatal(FARGS
, "Expanded ensemble not supported by MiMiC.");
176 gmx_fatal(FARGS
, "Simulated tempering not supported by MiMiC.");
180 gmx_fatal(FARGS
, "AWH not supported by MiMiC.");
182 if (replExParams
.exchangeInterval
> 0)
184 gmx_fatal(FARGS
, "Replica exchange not supported by MiMiC.");
186 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
188 gmx_fatal(FARGS
, "Essential dynamics not supported by MiMiC.");
192 gmx_fatal(FARGS
, "Interactive MD not supported by MiMiC.");
196 gmx_fatal(FARGS
, "Multiple simulations not supported by MiMiC.");
198 if (std::any_of(ir
->opts
.annealing
, ir
->opts
.annealing
+ ir
->opts
.ngtc
,
199 [](int i
) { return i
!= eannNO
; }))
201 gmx_fatal(FARGS
, "Simulated annealing not supported by MiMiC.");
204 /* Settings for rerun */
206 ir
->nstcalcenergy
= 1;
207 int nstglobalcomm
= 1;
208 const bool bNS
= true;
212 MimicCommunicator::init();
213 auto nonConstGlobalTopology
= const_cast<gmx_mtop_t
*>(top_global
);
214 MimicCommunicator::sendInitData(nonConstGlobalTopology
, state_global
->x
);
215 ir
->nsteps
= MimicCommunicator::getStepNumber();
218 ir
->nstxout_compressed
= 0;
219 const SimulationGroups
* groups
= &top_global
->groups
;
221 auto nonConstGlobalTopology
= const_cast<gmx_mtop_t
*>(top_global
);
222 nonConstGlobalTopology
->intermolecularExclusionGroup
= genQmmmIndices(*top_global
);
225 initialize_lambdas(fplog
, *ir
, MASTER(cr
), &state_global
->fep_state
, state_global
->lambda
, lam0
);
227 const bool simulationsShareState
= false;
228 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
229 mdModulesNotifier
, ir
, top_global
, oenv
, wcycle
,
230 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
231 gmx::EnergyOutput
energyOutput(mdoutf_get_fp_ene(outf
), top_global
, ir
, pull_work
,
232 mdoutf_get_fp_dhdl(outf
), true, StartingBehavior::NewSimulation
,
235 gstat
= global_stat_init(ir
);
237 /* Check for polarizable models and flexible constraints */
238 shellfc
= init_shell_flexcon(fplog
, top_global
, constr
? constr
->numFlexibleConstraints() : 0,
239 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
242 double io
= compute_io(ir
, top_global
->natoms
, *groups
, energyOutput
.numEnergyTerms(), 1);
243 if ((io
> 2000) && MASTER(cr
))
245 fprintf(stderr
, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io
);
249 // Local state only becomes valid now.
250 std::unique_ptr
<t_state
> stateInstance
;
253 gmx_localtop_t
top(top_global
->ffparams
);
255 if (DOMAINDECOMP(cr
))
257 stateInstance
= std::make_unique
<t_state
>();
258 state
= stateInstance
.get();
259 dd_init_local_state(cr
->dd
, state_global
, state
);
261 /* Distribute the charge groups over the nodes from the master node */
262 dd_partition_system(fplog
, mdlog
, ir
->init_step
, cr
, TRUE
, 1, state_global
, *top_global
, ir
,
263 imdSession
, pull_work
, state
, &f
, mdAtoms
, &top
, fr
, vsite
, constr
,
264 nrnb
, nullptr, FALSE
);
265 shouldCheckNumberOfBondedInteractions
= true;
266 gmx_bcast(sizeof(ir
->nsteps
), &ir
->nsteps
, cr
->mpi_comm_mygroup
);
270 state_change_natoms(state_global
, state_global
->natoms
);
271 /* Copy the pointer to the global state */
272 state
= state_global
;
274 mdAlgorithmsSetupAtomData(cr
, ir
, *top_global
, &top
, fr
, &f
, mdAtoms
, constr
, vsite
, shellfc
);
277 auto mdatoms
= mdAtoms
->mdatoms();
279 // NOTE: The global state is no longer used at this point.
280 // But state_global is still used as temporary storage space for writing
281 // the global state to file and potentially for replica exchange.
282 // (Global topology should persist.)
284 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
286 if (ir
->efep
!= efepNO
&& ir
->fepvals
->nstdhdl
!= 0)
288 doFreeEnergyPerturbation
= true;
294 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
295 bool bSumEkinhOld
= false;
296 t_vcm
* vcm
= nullptr;
297 compute_globals(gstat
, cr
, ir
, fr
, ekind
, makeConstArrayRef(state
->x
),
298 makeConstArrayRef(state
->v
), state
->box
, mdatoms
, nrnb
, vcm
, nullptr, enerd
,
299 force_vir
, shake_vir
, total_vir
, pres
, constr
, &nullSignaller
, state
->box
,
300 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
);
302 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
, top_global
, &top
,
303 makeConstArrayRef(state
->x
), state
->box
,
304 &shouldCheckNumberOfBondedInteractions
);
308 fprintf(stderr
, "starting MiMiC MD run '%s'\n\n", *(top_global
->name
));
309 if (mdrunOptions
.verbose
)
312 "Calculated time to finish depends on nsteps from "
313 "run input file,\nwhich may not correspond to the time "
314 "needed to process input trajectory.\n\n");
316 fprintf(fplog
, "\n");
319 walltime_accounting_start_time(walltime_accounting
);
320 wallcycle_start(wcycle
, ewcRUN
);
321 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
323 /***********************************************************
327 ************************************************************/
334 "Simulations has constraints. Constraints will "
335 "be handled by CPMD.");
341 "MiMiC does not report kinetic energy, total energy, temperature, virial and "
344 auto stopHandler
= stopHandlerBuilder
->getStopHandlerMD(
345 compat::not_null
<SimulationSignal
*>(&signals
[eglsSTOPCOND
]), false, MASTER(cr
),
346 ir
->nstlist
, mdrunOptions
.reproducible
, nstglobalcomm
, mdrunOptions
.maximumHoursToRun
,
347 ir
->nstlist
== 0, fplog
, step
, bNS
, walltime_accounting
);
349 // we don't do counter resetting in rerun - finish will always be valid
350 walltime_accounting_set_valid_finish(walltime_accounting
);
352 const DDBalanceRegionHandler
ddBalanceRegionHandler(cr
);
354 step
= ir
->init_step
;
357 /* and stop now if we should */
358 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
361 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
== ir
->nsteps
));
362 wallcycle_start(wcycle
, ewcSTEP
);
368 MimicCommunicator::getCoords(&state_global
->x
, state_global
->natoms
);
371 if (ir
->efep
!= efepNO
)
373 setCurrentLambdasLocal(step
, ir
->fepvals
, lam0
, state
->lambda
, state
->fep_state
);
378 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
379 if (constructVsites
&& DOMAINDECOMP(cr
))
382 "Vsite recalculation with -rerun is not implemented with domain "
384 "use a single rank");
388 if (DOMAINDECOMP(cr
))
390 /* Repartition the domain decomposition */
391 const bool bMasterState
= true;
392 dd_partition_system(fplog
, mdlog
, step
, cr
, bMasterState
, nstglobalcomm
, state_global
,
393 *top_global
, ir
, imdSession
, pull_work
, state
, &f
, mdAtoms
, &top
,
394 fr
, vsite
, constr
, nrnb
, wcycle
, mdrunOptions
.verbose
);
395 shouldCheckNumberOfBondedInteractions
= true;
400 EnergyOutput::printHeader(fplog
, step
, t
); /* can we improve the information printed here? */
403 if (ir
->efep
!= efepNO
)
405 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
408 force_flags
= (GMX_FORCE_STATECHANGED
| GMX_FORCE_DYNAMICBOX
| GMX_FORCE_ALLFORCES
409 | GMX_FORCE_VIRIAL
| // TODO: Get rid of this once #2649 is solved
410 GMX_FORCE_ENERGY
| (doFreeEnergyPerturbation
? GMX_FORCE_DHDL
: 0));
414 /* Now is the time to relax the shells */
415 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
, enforcedRotation
, step
, ir
,
416 imdSession
, pull_work
, bNS
, force_flags
, &top
, constr
, enerd
,
417 state
->natoms
, state
->x
.arrayRefWithPadding(),
418 state
->v
.arrayRefWithPadding(), state
->box
, state
->lambda
, &state
->hist
,
419 f
.arrayRefWithPadding(), force_vir
, mdatoms
, nrnb
, wcycle
, shellfc
,
420 fr
, runScheduleWork
, t
, mu_tot
, vsite
, ddBalanceRegionHandler
);
424 /* The coordinates (x) are shifted (to get whole molecules)
426 * This is parallellized as well, and does communication too.
427 * Check comments in sim_util.c
430 gmx_edsam
* ed
= nullptr;
431 do_force(fplog
, cr
, ms
, ir
, awh
, enforcedRotation
, imdSession
, pull_work
, step
, nrnb
,
432 wcycle
, &top
, state
->box
, state
->x
.arrayRefWithPadding(), &state
->hist
,
433 f
.arrayRefWithPadding(), force_vir
, mdatoms
, enerd
, state
->lambda
, fr
, runScheduleWork
,
434 vsite
, mu_tot
, t
, ed
, GMX_FORCE_NS
| force_flags
, ddBalanceRegionHandler
);
437 /* Now we have the energies and forces corresponding to the
438 * coordinates at time t.
441 const bool isCheckpointingStep
= false;
442 const bool doRerun
= false;
443 const bool bSumEkinhOld
= false;
444 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
, ir
, state
,
445 state_global
, observablesHistory
, top_global
, fr
, outf
,
446 energyOutput
, ekind
, f
, isCheckpointingStep
, doRerun
,
447 isLastStep
, mdrunOptions
.writeConfout
, bSumEkinhOld
);
450 stopHandler
->setSignal();
452 if (vsite
!= nullptr)
454 wallcycle_start(wcycle
, ewcVSITECONSTR
);
455 vsite
->construct(state
->x
, ir
->delta_t
, state
->v
, state
->box
);
456 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
460 const bool doInterSimSignal
= false;
461 const bool doIntraSimSignal
= true;
462 bool bSumEkinhOld
= false;
463 t_vcm
* vcm
= nullptr;
464 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
466 compute_globals(gstat
, cr
, ir
, fr
, ekind
, makeConstArrayRef(state
->x
),
467 makeConstArrayRef(state
->v
), state
->box
, mdatoms
, nrnb
, vcm
, wcycle
,
468 enerd
, nullptr, nullptr, nullptr, nullptr, constr
, &signaller
,
469 state
->box
, &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
470 CGLO_GSTAT
| CGLO_ENERGY
471 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
473 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
, top_global
,
474 &top
, makeConstArrayRef(state
->x
), state
->box
,
475 &shouldCheckNumberOfBondedInteractions
);
479 gmx::HostVector
<gmx::RVec
> fglobal(top_global
->natoms
);
480 gmx::ArrayRef
<gmx::RVec
> ftemp
;
481 gmx::ArrayRef
<const gmx::RVec
> flocal
= gmx::makeArrayRef(f
);
482 if (DOMAINDECOMP(cr
))
484 ftemp
= gmx::makeArrayRef(fglobal
);
485 dd_collect_vec(cr
->dd
, state
, flocal
, ftemp
);
489 ftemp
= gmx::makeArrayRef(f
);
494 MimicCommunicator::sendEnergies(enerd
->term
[F_EPOT
]);
495 MimicCommunicator::sendForces(ftemp
, state_global
->natoms
);
500 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
501 the virial that should probably be addressed eventually. state->veta has better properies,
502 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
503 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
505 if (ir
->efep
!= efepNO
)
507 /* Sum up the foreign energy and dhdl terms for md and sd.
508 Currently done every step so that dhdl is correct in the .edr */
509 accumulateKineticLambdaComponents(enerd
, state
->lambda
, *ir
->fepvals
);
515 const bool bCalcEnerStep
= true;
516 energyOutput
.addDataAtEnergyStep(doFreeEnergyPerturbation
, bCalcEnerStep
, t
,
517 mdatoms
->tmass
, enerd
, state
, ir
->fepvals
,
518 ir
->expandedvals
, state
->box
, shake_vir
, force_vir
,
519 total_vir
, pres
, ekind
, mu_tot
, constr
);
521 const bool do_ene
= true;
522 const bool do_log
= true;
524 const bool do_dr
= ir
->nstdisreout
!= 0;
525 const bool do_or
= ir
->nstorireout
!= 0;
527 EnergyOutput::printAnnealingTemperatures(do_log
? fplog
: nullptr, groups
, &(ir
->opts
));
528 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
,
529 do_log
? fplog
: nullptr, step
, t
,
530 &fr
->listedForces
->fcdata(), awh
);
532 if (do_per_step(step
, ir
->nstlog
))
534 if (fflush(fplog
) != 0)
536 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
541 /* Print the remaining wall clock time for the run */
542 if (isMasterSimMasterRank(ms
, MASTER(cr
)) && (mdrunOptions
.verbose
|| gmx_got_usr_signal()))
546 fprintf(stderr
, "\n");
548 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
551 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
552 if (DOMAINDECOMP(cr
) && wcycle
)
554 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
557 /* increase the MD step number */
561 /* End of main MD loop */
563 /* Closing TNG files can include compressing data. Therefore it is good to do that
564 * before stopping the time measurements. */
565 mdoutf_tng_close(outf
);
567 /* Stop measuring walltime */
568 walltime_accounting_end_time(walltime_accounting
);
572 MimicCommunicator::finalize();
575 if (!thisRankHasDuty(cr
, DUTY_PME
))
577 /* Tell the PME only node to finish */
578 gmx_pme_send_finish(cr
);
583 done_shellfc(fplog
, shellfc
, step_rel
);
585 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);