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/applied_forces/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/freeenergyparameters.h"
85 #include "gromacs/mdlib/md_support.h"
86 #include "gromacs/mdlib/mdatoms.h"
87 #include "gromacs/mdlib/mdoutf.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/resethandler.h"
90 #include "gromacs/mdlib/sighandler.h"
91 #include "gromacs/mdlib/simulationsignal.h"
92 #include "gromacs/mdlib/stat.h"
93 #include "gromacs/mdlib/stophandler.h"
94 #include "gromacs/mdlib/tgroup.h"
95 #include "gromacs/mdlib/trajectory_writing.h"
96 #include "gromacs/mdlib/update.h"
97 #include "gromacs/mdlib/vcm.h"
98 #include "gromacs/mdlib/vsite.h"
99 #include "gromacs/mdrunutility/handlerestart.h"
100 #include "gromacs/mdrunutility/multisim.h"
101 #include "gromacs/mdrunutility/printtime.h"
102 #include "gromacs/mdtypes/awh_history.h"
103 #include "gromacs/mdtypes/awh_params.h"
104 #include "gromacs/mdtypes/commrec.h"
105 #include "gromacs/mdtypes/df_history.h"
106 #include "gromacs/mdtypes/enerdata.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/forcebuffers.h"
109 #include "gromacs/mdtypes/forcerec.h"
110 #include "gromacs/mdtypes/group.h"
111 #include "gromacs/mdtypes/inputrec.h"
112 #include "gromacs/mdtypes/interaction_const.h"
113 #include "gromacs/mdtypes/md_enums.h"
114 #include "gromacs/mdtypes/mdatom.h"
115 #include "gromacs/mdtypes/mdrunoptions.h"
116 #include "gromacs/mdtypes/observableshistory.h"
117 #include "gromacs/mdtypes/simulation_workload.h"
118 #include "gromacs/mdtypes/state.h"
119 #include "gromacs/mimic/communicator.h"
120 #include "gromacs/mimic/utilities.h"
121 #include "gromacs/pbcutil/pbc.h"
122 #include "gromacs/pulling/pull.h"
123 #include "gromacs/timing/wallcycle.h"
124 #include "gromacs/timing/walltime_accounting.h"
125 #include "gromacs/topology/atoms.h"
126 #include "gromacs/topology/idef.h"
127 #include "gromacs/topology/mtop_util.h"
128 #include "gromacs/topology/topology.h"
129 #include "gromacs/trajectory/trajectoryframe.h"
130 #include "gromacs/utility/basedefinitions.h"
131 #include "gromacs/utility/cstringutil.h"
132 #include "gromacs/utility/fatalerror.h"
133 #include "gromacs/utility/logger.h"
134 #include "gromacs/utility/real.h"
136 #include "legacysimulator.h"
137 #include "replicaexchange.h"
140 using gmx::SimulationSignaller
;
142 void gmx::LegacySimulator::do_mimic()
144 t_inputrec
* ir
= inputrec
;
145 int64_t step
, step_rel
;
147 bool isLastStep
= false;
148 bool doFreeEnergyPerturbation
= false;
149 unsigned int force_flags
;
150 tensor force_vir
, shake_vir
, total_vir
, pres
;
153 gmx_global_stat_t gstat
;
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;
214 MimicCommunicator::init();
215 auto nonConstGlobalTopology
= const_cast<gmx_mtop_t
*>(top_global
);
216 MimicCommunicator::sendInitData(nonConstGlobalTopology
, state_global
->x
);
217 ir
->nsteps
= MimicCommunicator::getStepNumber();
220 ir
->nstxout_compressed
= 0;
221 const SimulationGroups
* groups
= &top_global
->groups
;
223 auto nonConstGlobalTopology
= const_cast<gmx_mtop_t
*>(top_global
);
224 nonConstGlobalTopology
->intermolecularExclusionGroup
= genQmmmIndices(*top_global
);
227 initialize_lambdas(fplog
, *ir
, MASTER(cr
), &state_global
->fep_state
, state_global
->lambda
);
229 const bool simulationsShareState
= false;
230 gmx_mdoutf
* outf
= init_mdoutf(fplog
, nfile
, fnm
, mdrunOptions
, cr
, outputProvider
,
231 mdModulesNotifier
, ir
, top_global
, oenv
, wcycle
,
232 StartingBehavior::NewSimulation
, simulationsShareState
, ms
);
233 gmx::EnergyOutput
energyOutput(mdoutf_get_fp_ene(outf
), top_global
, ir
, pull_work
,
234 mdoutf_get_fp_dhdl(outf
), true, StartingBehavior::NewSimulation
,
235 simulationsShareState
, mdModulesNotifier
);
237 gstat
= global_stat_init(ir
);
239 /* Check for polarizable models and flexible constraints */
240 shellfc
= init_shell_flexcon(fplog
, top_global
, constr
? constr
->numFlexibleConstraints() : 0,
241 ir
->nstcalcenergy
, DOMAINDECOMP(cr
),
242 runScheduleWork
->simulationWork
.useGpuPme
);
245 double io
= compute_io(ir
, top_global
->natoms
, *groups
, energyOutput
.numEnergyTerms(), 1);
246 if ((io
> 2000) && MASTER(cr
))
248 fprintf(stderr
, "\nWARNING: This run will generate roughly %.0f Mb of data\n\n", io
);
252 // Local state only becomes valid now.
253 std::unique_ptr
<t_state
> stateInstance
;
256 gmx_localtop_t
top(top_global
->ffparams
);
258 if (DOMAINDECOMP(cr
))
260 stateInstance
= std::make_unique
<t_state
>();
261 state
= stateInstance
.get();
262 dd_init_local_state(cr
->dd
, state_global
, state
);
264 /* Distribute the charge groups over the nodes from the master node */
265 dd_partition_system(fplog
, mdlog
, ir
->init_step
, cr
, TRUE
, 1, state_global
, *top_global
, ir
,
266 imdSession
, pull_work
, state
, &f
, mdAtoms
, &top
, fr
, vsite
, constr
,
267 nrnb
, nullptr, FALSE
);
268 shouldCheckNumberOfBondedInteractions
= true;
269 gmx_bcast(sizeof(ir
->nsteps
), &ir
->nsteps
, cr
->mpi_comm_mygroup
);
273 state_change_natoms(state_global
, state_global
->natoms
);
274 /* Copy the pointer to the global state */
275 state
= state_global
;
277 mdAlgorithmsSetupAtomData(cr
, ir
, *top_global
, &top
, fr
, &f
, mdAtoms
, constr
, vsite
, shellfc
);
280 auto mdatoms
= mdAtoms
->mdatoms();
282 // NOTE: The global state is no longer used at this point.
283 // But state_global is still used as temporary storage space for writing
284 // the global state to file and potentially for replica exchange.
285 // (Global topology should persist.)
287 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
289 if (ir
->efep
!= efepNO
&& ir
->fepvals
->nstdhdl
!= 0)
291 doFreeEnergyPerturbation
= true;
297 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
298 bool bSumEkinhOld
= false;
299 t_vcm
* vcm
= nullptr;
300 compute_globals(gstat
, cr
, ir
, fr
, ekind
, makeConstArrayRef(state
->x
),
301 makeConstArrayRef(state
->v
), state
->box
, mdatoms
, nrnb
, vcm
, nullptr, enerd
,
302 force_vir
, shake_vir
, total_vir
, pres
, constr
, &nullSignaller
, state
->box
,
303 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
);
305 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
, top_global
, &top
,
306 makeConstArrayRef(state
->x
), state
->box
,
307 &shouldCheckNumberOfBondedInteractions
);
311 fprintf(stderr
, "starting MiMiC MD run '%s'\n\n", *(top_global
->name
));
312 if (mdrunOptions
.verbose
)
315 "Calculated time to finish depends on nsteps from "
316 "run input file,\nwhich may not correspond to the time "
317 "needed to process input trajectory.\n\n");
319 fprintf(fplog
, "\n");
322 walltime_accounting_start_time(walltime_accounting
);
323 wallcycle_start(wcycle
, ewcRUN
);
324 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
326 /***********************************************************
330 ************************************************************/
337 "Simulations has constraints. Constraints will "
338 "be handled by CPMD.");
344 "MiMiC does not report kinetic energy, total energy, temperature, virial and "
347 step
= ir
->init_step
;
350 auto stopHandler
= stopHandlerBuilder
->getStopHandlerMD(
351 compat::not_null
<SimulationSignal
*>(&signals
[eglsSTOPCOND
]), false, MASTER(cr
),
352 ir
->nstlist
, mdrunOptions
.reproducible
, nstglobalcomm
, mdrunOptions
.maximumHoursToRun
,
353 ir
->nstlist
== 0, fplog
, step
, bNS
, walltime_accounting
);
355 // we don't do counter resetting in rerun - finish will always be valid
356 walltime_accounting_set_valid_finish(walltime_accounting
);
358 const DDBalanceRegionHandler
ddBalanceRegionHandler(cr
);
360 /* and stop now if we should */
361 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
364 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
== ir
->nsteps
));
365 wallcycle_start(wcycle
, ewcSTEP
);
371 MimicCommunicator::getCoords(&state_global
->x
, state_global
->natoms
);
374 if (ir
->efep
!= efepNO
)
376 state
->lambda
= currentLambdas(step
, *(ir
->fepvals
), state_global
->fep_state
);
381 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
382 if (constructVsites
&& DOMAINDECOMP(cr
))
385 "Vsite recalculation with -rerun is not implemented with domain "
387 "use a single rank");
391 if (DOMAINDECOMP(cr
))
393 /* Repartition the domain decomposition */
394 const bool bMasterState
= true;
395 dd_partition_system(fplog
, mdlog
, step
, cr
, bMasterState
, nstglobalcomm
, state_global
,
396 *top_global
, ir
, imdSession
, pull_work
, state
, &f
, mdAtoms
, &top
,
397 fr
, vsite
, constr
, nrnb
, wcycle
, mdrunOptions
.verbose
);
398 shouldCheckNumberOfBondedInteractions
= true;
403 EnergyOutput::printHeader(fplog
, step
, t
); /* can we improve the information printed here? */
406 if (ir
->efep
!= efepNO
)
408 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
411 force_flags
= (GMX_FORCE_STATECHANGED
| GMX_FORCE_DYNAMICBOX
| GMX_FORCE_ALLFORCES
412 | GMX_FORCE_VIRIAL
| // TODO: Get rid of this once #2649 is solved
413 GMX_FORCE_ENERGY
| (doFreeEnergyPerturbation
? GMX_FORCE_DHDL
: 0));
417 /* Now is the time to relax the shells */
418 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
, enforcedRotation
, step
, ir
,
419 imdSession
, pull_work
, bNS
, force_flags
, &top
, constr
, enerd
,
420 state
->natoms
, state
->x
.arrayRefWithPadding(),
421 state
->v
.arrayRefWithPadding(), state
->box
, state
->lambda
,
422 &state
->hist
, &f
.view(), force_vir
, mdatoms
, nrnb
, wcycle
, shellfc
,
423 fr
, runScheduleWork
, t
, mu_tot
, vsite
, ddBalanceRegionHandler
);
427 /* The coordinates (x) are shifted (to get whole molecules)
429 * This is parallellized as well, and does communication too.
430 * Check comments in sim_util.c
433 gmx_edsam
* ed
= nullptr;
434 do_force(fplog
, cr
, ms
, ir
, awh
, enforcedRotation
, imdSession
, pull_work
, step
, nrnb
,
435 wcycle
, &top
, state
->box
, state
->x
.arrayRefWithPadding(), &state
->hist
,
436 &f
.view(), force_vir
, mdatoms
, enerd
, state
->lambda
, fr
, runScheduleWork
,
437 vsite
, mu_tot
, t
, ed
, GMX_FORCE_NS
| force_flags
, ddBalanceRegionHandler
);
440 /* Now we have the energies and forces corresponding to the
441 * coordinates at time t.
444 const bool isCheckpointingStep
= false;
445 const bool doRerun
= false;
446 const bool bSumEkinhOld
= false;
447 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
, ir
, state
,
448 state_global
, observablesHistory
, top_global
, fr
, outf
,
449 energyOutput
, ekind
, f
.view().force(), isCheckpointingStep
,
450 doRerun
, isLastStep
, mdrunOptions
.writeConfout
, bSumEkinhOld
);
453 stopHandler
->setSignal();
455 if (vsite
!= nullptr)
457 wallcycle_start(wcycle
, ewcVSITECONSTR
);
458 vsite
->construct(state
->x
, ir
->delta_t
, state
->v
, state
->box
);
459 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
463 const bool doInterSimSignal
= false;
464 const bool doIntraSimSignal
= true;
465 bool bSumEkinhOld
= false;
466 t_vcm
* vcm
= nullptr;
467 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
469 compute_globals(gstat
, cr
, ir
, fr
, ekind
, makeConstArrayRef(state
->x
),
470 makeConstArrayRef(state
->v
), state
->box
, mdatoms
, nrnb
, vcm
, wcycle
,
471 enerd
, nullptr, nullptr, nullptr, nullptr, constr
, &signaller
,
472 state
->box
, &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
473 CGLO_GSTAT
| CGLO_ENERGY
474 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
476 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
, top_global
,
477 &top
, makeConstArrayRef(state
->x
), state
->box
,
478 &shouldCheckNumberOfBondedInteractions
);
482 gmx::HostVector
<gmx::RVec
> fglobal(top_global
->natoms
);
483 gmx::ArrayRef
<gmx::RVec
> ftemp
;
484 gmx::ArrayRef
<const gmx::RVec
> flocal
= f
.view().force();
485 if (DOMAINDECOMP(cr
))
487 ftemp
= gmx::makeArrayRef(fglobal
);
488 dd_collect_vec(cr
->dd
, state
->ddp_count
, state
->ddp_count_cg_gl
, state
->cg_gl
,
493 ftemp
= f
.view().force();
498 MimicCommunicator::sendEnergies(enerd
->term
[F_EPOT
]);
499 MimicCommunicator::sendForces(ftemp
, state_global
->natoms
);
504 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
505 the virial that should probably be addressed eventually. state->veta has better properies,
506 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
507 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
509 if (ir
->efep
!= efepNO
)
511 /* Sum up the foreign energy and dhdl terms for md and sd.
512 Currently done every step so that dhdl is correct in the .edr */
513 accumulateKineticLambdaComponents(enerd
, state
->lambda
, *ir
->fepvals
);
519 const bool bCalcEnerStep
= true;
520 energyOutput
.addDataAtEnergyStep(
521 doFreeEnergyPerturbation
, bCalcEnerStep
, t
, mdatoms
->tmass
, enerd
, ir
->fepvals
,
522 ir
->expandedvals
, state
->box
,
523 PTCouplingArrays({ state
->boxv
, state
->nosehoover_xi
, state
->nosehoover_vxi
,
524 state
->nhpres_xi
, state
->nhpres_vxi
}),
525 state
->fep_state
, shake_vir
, force_vir
, total_vir
, pres
, ekind
, mu_tot
, constr
);
527 const bool do_ene
= true;
528 const bool do_log
= true;
530 const bool do_dr
= ir
->nstdisreout
!= 0;
531 const bool do_or
= ir
->nstorireout
!= 0;
533 EnergyOutput::printAnnealingTemperatures(do_log
? fplog
: nullptr, groups
, &(ir
->opts
));
534 energyOutput
.printStepToEnergyFile(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
,
535 do_log
? fplog
: nullptr, step
, t
, fr
->fcdata
.get(), awh
);
537 if (do_per_step(step
, ir
->nstlog
))
539 if (fflush(fplog
) != 0)
541 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
546 /* Print the remaining wall clock time for the run */
547 if (isMasterSimMasterRank(ms
, MASTER(cr
)) && (mdrunOptions
.verbose
|| gmx_got_usr_signal()))
551 fprintf(stderr
, "\n");
553 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
556 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
557 if (DOMAINDECOMP(cr
) && wcycle
)
559 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
562 /* increase the MD step number */
566 /* End of main MD loop */
568 /* Closing TNG files can include compressing data. Therefore it is good to do that
569 * before stopping the time measurements. */
570 mdoutf_tng_close(outf
);
572 /* Stop measuring walltime */
573 walltime_accounting_end_time(walltime_accounting
);
577 MimicCommunicator::finalize();
580 if (!thisRankHasDuty(cr
, DUTY_PME
))
582 /* Tell the PME only node to finish */
583 gmx_pme_send_finish(cr
);
588 done_shellfc(fplog
, shellfc
, step_rel
);
590 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);