2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018, 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/compat/make_unique.h"
56 #include "gromacs/domdec/collect.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/partition.h"
61 #include "gromacs/essentialdynamics/edsam.h"
62 #include "gromacs/ewald/pme.h"
63 #include "gromacs/ewald/pme-load-balancing.h"
64 #include "gromacs/fileio/trxio.h"
65 #include "gromacs/gmxlib/network.h"
66 #include "gromacs/gmxlib/nrnb.h"
67 #include "gromacs/gpu_utils/gpu_utils.h"
68 #include "gromacs/imd/imd.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/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/mdebin.h"
85 #include "gromacs/mdlib/mdoutf.h"
86 #include "gromacs/mdlib/mdrun.h"
87 #include "gromacs/mdlib/mdsetup.h"
88 #include "gromacs/mdlib/membed.h"
89 #include "gromacs/mdlib/nb_verlet.h"
90 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
91 #include "gromacs/mdlib/ns.h"
92 #include "gromacs/mdlib/resethandler.h"
93 #include "gromacs/mdlib/shellfc.h"
94 #include "gromacs/mdlib/sighandler.h"
95 #include "gromacs/mdlib/sim_util.h"
96 #include "gromacs/mdlib/simulationsignal.h"
97 #include "gromacs/mdlib/stophandler.h"
98 #include "gromacs/mdlib/tgroup.h"
99 #include "gromacs/mdlib/trajectory_writing.h"
100 #include "gromacs/mdlib/update.h"
101 #include "gromacs/mdlib/vcm.h"
102 #include "gromacs/mdlib/vsite.h"
103 #include "gromacs/mdtypes/awh-history.h"
104 #include "gromacs/mdtypes/awh-params.h"
105 #include "gromacs/mdtypes/commrec.h"
106 #include "gromacs/mdtypes/df_history.h"
107 #include "gromacs/mdtypes/energyhistory.h"
108 #include "gromacs/mdtypes/fcdata.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/observableshistory.h"
116 #include "gromacs/mdtypes/state.h"
117 #include "gromacs/mimic/MimicCommunicator.h"
118 #include "gromacs/mimic/MimicUtils.h"
119 #include "gromacs/pbcutil/mshift.h"
120 #include "gromacs/pbcutil/pbc.h"
121 #include "gromacs/pulling/pull.h"
122 #include "gromacs/swap/swapcoords.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"
135 #include "gromacs/utility/smalloc.h"
137 #include "integrator.h"
138 #include "replicaexchange.h"
140 using gmx::SimulationSignaller
;
142 void gmx::Integrator::do_mimic()
144 t_inputrec
*ir
= inputrec
;
145 gmx_mdoutf
*outf
= nullptr;
146 int64_t step
, step_rel
;
147 double t
, lam0
[efptNR
];
148 bool isLastStep
= false;
149 bool doFreeEnergyPerturbation
= false;
151 tensor force_vir
, shake_vir
, total_vir
, pres
;
154 t_mdebin
*mdebin
= nullptr;
155 gmx_enerdata_t
*enerd
;
156 PaddedVector
<gmx::RVec
> f
{};
157 gmx_global_stat_t gstat
;
158 t_graph
*graph
= nullptr;
159 gmx_groups_t
*groups
;
160 gmx_ekindata_t
*ekind
;
161 gmx_shellfc_t
*shellfc
;
165 /* Domain decomposition could incorrectly miss a bonded
166 interaction, but checking for that requires a global
167 communication stage, which does not otherwise happen in DD
168 code. So we do that alongside the first global energy reduction
169 after a new DD is made. These variables handle whether the
170 check happens, and the result it returns. */
171 bool shouldCheckNumberOfBondedInteractions
= false;
172 int totalNumberOfBondedInteractions
= -1;
174 SimulationSignals signals
;
175 // Most global communnication stages don't propagate mdrun
176 // signals, and will use this object to achieve that.
177 SimulationSignaller
nullSignaller(nullptr, nullptr, nullptr, false, false);
181 gmx_fatal(FARGS
, "Expanded ensemble not supported by MiMiC.");
185 gmx_fatal(FARGS
, "Simulated tempering not supported by MiMiC.");
189 gmx_fatal(FARGS
, "AWH not supported by MiMiC.");
191 if (replExParams
.exchangeInterval
> 0)
193 gmx_fatal(FARGS
, "Replica exchange not supported by MiMiC.");
195 if (opt2bSet("-ei", nfile
, fnm
) || observablesHistory
->edsamHistory
!= nullptr)
197 gmx_fatal(FARGS
, "Essential dynamics not supported by MiMiC.");
201 gmx_fatal(FARGS
, "Interactive MD not supported by MiMiC.");
205 gmx_fatal(FARGS
, "Multiple simulations not supported by MiMiC.");
207 if (std::any_of(ir
->opts
.annealing
, ir
->opts
.annealing
+ ir
->opts
.ngtc
,
208 [](int i
){return i
!= eannNO
; }))
210 gmx_fatal(FARGS
, "Simulated annealing not supported by MiMiC.");
213 /* Settings for rerun */
215 ir
->nstcalcenergy
= 1;
216 int nstglobalcomm
= 1;
217 const bool bNS
= true;
219 // Communicator to interact with MiMiC
220 MimicCommunicator mimicCommunicator
{};
223 mimicCommunicator
.init();
224 mimicCommunicator
.sendInitData(top_global
, state_global
->x
);
225 ir
->nsteps
= mimicCommunicator
.getStepNumber();
228 ir
->nstxout_compressed
= 0;
229 groups
= &top_global
->groups
;
230 top_global
->intermolecularExclusionGroup
= genQmmmIndices(*top_global
);
233 init_rerun(fplog
, cr
, outputProvider
, ir
, oenv
, mdrunOptions
,
234 state_global
, lam0
, nrnb
, top_global
,
235 nfile
, fnm
, &outf
, &mdebin
, wcycle
);
237 /* Energy terms and groups */
239 init_enerdata(top_global
->groups
.grps
[egcENER
].nr
, ir
->fepvals
->n_lambda
,
242 /* Kinetic energy data */
244 init_ekindata(fplog
, top_global
, &(ir
->opts
), ekind
);
245 /* Copy the cos acceleration to the groups struct */
246 ekind
->cosacc
.cos_accel
= ir
->cos_accel
;
248 gstat
= global_stat_init(ir
);
250 /* Check for polarizable models and flexible constraints */
251 shellfc
= init_shell_flexcon(fplog
,
252 top_global
, constr
? constr
->numFlexibleConstraints() : 0,
253 ir
->nstcalcenergy
, DOMAINDECOMP(cr
));
256 double io
= compute_io(ir
, top_global
->natoms
, groups
, mdebin
->ebin
->nener
, 1);
257 if ((io
> 2000) && MASTER(cr
))
260 "\nWARNING: This run will generate roughly %.0f Mb of data\n\n",
265 // Local state only becomes valid now.
266 std::unique_ptr
<t_state
> stateInstance
;
269 if (DOMAINDECOMP(cr
))
271 top
= dd_init_local_top(top_global
);
273 stateInstance
= compat::make_unique
<t_state
>();
274 state
= stateInstance
.get();
275 dd_init_local_state(cr
->dd
, state_global
, state
);
277 /* Distribute the charge groups over the nodes from the master node */
278 dd_partition_system(fplog
, mdlog
, ir
->init_step
, cr
, TRUE
, 1,
279 state_global
, top_global
, ir
,
280 state
, &f
, mdAtoms
, top
, fr
,
282 nrnb
, nullptr, FALSE
);
283 shouldCheckNumberOfBondedInteractions
= true;
284 gmx_bcast(sizeof(ir
->nsteps
), &ir
->nsteps
, cr
);
288 state_change_natoms(state_global
, state_global
->natoms
);
289 /* We need to allocate one element extra, since we might use
290 * (unaligned) 4-wide SIMD loads to access rvec entries.
292 f
.resizeWithPadding(state_global
->natoms
);
293 /* Copy the pointer to the global state */
294 state
= state_global
;
297 mdAlgorithmsSetupAtomData(cr
, ir
, top_global
, top
, fr
,
298 &graph
, mdAtoms
, constr
, vsite
, shellfc
);
301 auto mdatoms
= mdAtoms
->mdatoms();
303 // NOTE: The global state is no longer used at this point.
304 // But state_global is still used as temporary storage space for writing
305 // the global state to file and potentially for replica exchange.
306 // (Global topology should persist.)
308 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
310 if (ir
->efep
!= efepNO
&& ir
->fepvals
->nstdhdl
!= 0)
312 doFreeEnergyPerturbation
= true;
316 int cglo_flags
= (CGLO_INITIALIZATION
| CGLO_GSTAT
|
317 (shouldCheckNumberOfBondedInteractions
?
318 CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0));
319 bool bSumEkinhOld
= false;
320 t_vcm
*vcm
= nullptr;
321 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
322 nullptr, enerd
, force_vir
, shake_vir
, total_vir
, pres
, mu_tot
,
323 constr
, &nullSignaller
, state
->box
,
324 &totalNumberOfBondedInteractions
, &bSumEkinhOld
, cglo_flags
);
326 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
327 top_global
, top
, state
,
328 &shouldCheckNumberOfBondedInteractions
);
332 fprintf(stderr
, "starting MiMiC MD run '%s'\n\n",
333 *(top_global
->name
));
334 if (mdrunOptions
.verbose
)
336 fprintf(stderr
, "Calculated time to finish depends on nsteps from "
337 "run input file,\nwhich may not correspond to the time "
338 "needed to process input trajectory.\n\n");
340 fprintf(fplog
, "\n");
343 walltime_accounting_start_time(walltime_accounting
);
344 wallcycle_start(wcycle
, ewcRUN
);
345 print_start(fplog
, cr
, walltime_accounting
, "mdrun");
347 /***********************************************************
351 ************************************************************/
355 GMX_LOG(mdlog
.info
).asParagraph().
356 appendText("Simulations has constraints. Constraints will "
357 "be handled by CPMD.");
360 GMX_LOG(mdlog
.info
).asParagraph().
361 appendText("MiMiC does not report kinetic energy, total energy, temperature, virial and pressure.");
363 auto stopHandler
= stopHandlerBuilder
->getStopHandlerMD(
364 compat::not_null
<SimulationSignal
*>(&signals
[eglsSTOPCOND
]), false,
365 MASTER(cr
), ir
->nstlist
, mdrunOptions
.reproducible
, nstglobalcomm
,
366 mdrunOptions
.maximumHoursToRun
, ir
->nstlist
== 0, fplog
, step
, bNS
, walltime_accounting
);
368 // we don't do counter resetting in rerun - finish will always be valid
369 walltime_accounting_set_valid_finish(walltime_accounting
);
371 DdOpenBalanceRegionBeforeForceComputation ddOpenBalanceRegion
=
373 DdOpenBalanceRegionBeforeForceComputation::yes
:
374 DdOpenBalanceRegionBeforeForceComputation::no
);
375 DdCloseBalanceRegionAfterForceComputation ddCloseBalanceRegion
=
377 DdCloseBalanceRegionAfterForceComputation::yes
:
378 DdCloseBalanceRegionAfterForceComputation::no
);
380 step
= ir
->init_step
;
383 /* and stop now if we should */
384 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
> ir
->nsteps
));
387 isLastStep
= (isLastStep
|| (ir
->nsteps
>= 0 && step_rel
== ir
->nsteps
));
388 wallcycle_start(wcycle
, ewcSTEP
);
394 mimicCommunicator
.getCoords(&state_global
->x
, state_global
->natoms
);
397 if (ir
->efep
!= efepNO
)
399 setCurrentLambdasLocal(step
, ir
->fepvals
, lam0
, state
);
404 const bool constructVsites
= ((vsite
!= nullptr) && mdrunOptions
.rerunConstructVsites
);
405 if (constructVsites
&& DOMAINDECOMP(cr
))
407 gmx_fatal(FARGS
, "Vsite recalculation with -rerun is not implemented with domain decomposition, "
408 "use a single rank");
412 if (DOMAINDECOMP(cr
))
414 /* Repartition the domain decomposition */
415 const bool bMasterState
= true;
416 dd_partition_system(fplog
, mdlog
, step
, cr
,
417 bMasterState
, nstglobalcomm
,
418 state_global
, top_global
, ir
,
419 state
, &f
, mdAtoms
, top
, fr
,
422 mdrunOptions
.verbose
);
423 shouldCheckNumberOfBondedInteractions
= true;
428 print_ebin_header(fplog
, step
, t
); /* can we improve the information printed here? */
431 if (ir
->efep
!= efepNO
)
433 update_mdatoms(mdatoms
, state
->lambda
[efptMASS
]);
436 force_flags
= (GMX_FORCE_STATECHANGED
|
437 GMX_FORCE_DYNAMICBOX
|
438 GMX_FORCE_ALLFORCES
|
439 GMX_FORCE_VIRIAL
| // TODO: Get rid of this once #2649 is solved
441 (doFreeEnergyPerturbation
? GMX_FORCE_DHDL
: 0));
445 /* Now is the time to relax the shells */
446 relax_shell_flexcon(fplog
, cr
, ms
, mdrunOptions
.verbose
,
447 enforcedRotation
, step
,
448 ir
, bNS
, force_flags
, top
,
450 state
, f
.arrayRefWithPadding(), force_vir
, mdatoms
,
451 nrnb
, wcycle
, graph
, groups
,
452 shellfc
, fr
, t
, mu_tot
,
454 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
458 /* The coordinates (x) are shifted (to get whole molecules)
460 * This is parallellized as well, and does communication too.
461 * Check comments in sim_util.c
464 gmx_edsam
*ed
= nullptr;
465 do_force(fplog
, cr
, ms
, ir
, awh
, enforcedRotation
,
466 step
, nrnb
, wcycle
, top
, groups
,
467 state
->box
, state
->x
.arrayRefWithPadding(), &state
->hist
,
468 f
.arrayRefWithPadding(), force_vir
, mdatoms
, enerd
, fcd
,
469 state
->lambda
, graph
,
470 fr
, vsite
, mu_tot
, t
, ed
,
471 GMX_FORCE_NS
| force_flags
,
472 ddOpenBalanceRegion
, ddCloseBalanceRegion
);
475 /* Now we have the energies and forces corresponding to the
476 * coordinates at time t.
479 const bool isCheckpointingStep
= false;
480 const bool doRerun
= false;
481 const bool bSumEkinhOld
= false;
482 do_md_trajectory_writing(fplog
, cr
, nfile
, fnm
, step
, step_rel
, t
,
483 ir
, state
, state_global
, observablesHistory
,
485 outf
, mdebin
, ekind
, f
,
486 isCheckpointingStep
, doRerun
, isLastStep
,
487 mdrunOptions
.writeConfout
,
491 stopHandler
->setSignal();
495 /* Need to unshift here */
496 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
499 if (vsite
!= nullptr)
501 wallcycle_start(wcycle
, ewcVSITECONSTR
);
502 if (graph
!= nullptr)
504 shift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
506 construct_vsites(vsite
, as_rvec_array(state
->x
.data()), ir
->delta_t
, as_rvec_array(state
->v
.data()),
507 top
->idef
.iparams
, top
->idef
.il
,
508 fr
->ePBC
, fr
->bMolPBC
, cr
, state
->box
);
510 if (graph
!= nullptr)
512 unshift_self(graph
, state
->box
, as_rvec_array(state
->x
.data()));
514 wallcycle_stop(wcycle
, ewcVSITECONSTR
);
518 const bool doInterSimSignal
= false;
519 const bool doIntraSimSignal
= true;
520 bool bSumEkinhOld
= false;
521 t_vcm
*vcm
= nullptr;
522 SimulationSignaller
signaller(&signals
, cr
, ms
, doInterSimSignal
, doIntraSimSignal
);
524 compute_globals(fplog
, gstat
, cr
, ir
, fr
, ekind
, state
, mdatoms
, nrnb
, vcm
,
525 wcycle
, enerd
, nullptr, nullptr, nullptr, nullptr, mu_tot
,
528 &totalNumberOfBondedInteractions
, &bSumEkinhOld
,
529 CGLO_GSTAT
| CGLO_ENERGY
530 | (shouldCheckNumberOfBondedInteractions
? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS
: 0)
532 checkNumberOfBondedInteractions(mdlog
, cr
, totalNumberOfBondedInteractions
,
533 top_global
, top
, state
,
534 &shouldCheckNumberOfBondedInteractions
);
538 gmx::HostVector
<gmx::RVec
> fglobal(top_global
->natoms
);
539 gmx::ArrayRef
<gmx::RVec
> ftemp
;
540 gmx::ArrayRef
<const gmx::RVec
> flocal
= gmx::makeArrayRef(f
);
541 if (DOMAINDECOMP(cr
))
543 ftemp
= gmx::makeArrayRef(fglobal
);
544 dd_collect_vec(cr
->dd
, state
, flocal
, ftemp
);
548 ftemp
= gmx::makeArrayRef(f
);
553 mimicCommunicator
.sendEnergies(enerd
->term
[F_EPOT
]);
554 mimicCommunicator
.sendForces(ftemp
, state_global
->natoms
);
560 /* Note: this is OK, but there are some numerical precision issues with using the convergence of
561 the virial that should probably be addressed eventually. state->veta has better properies,
562 but what we actually need entering the new cycle is the new shake_vir value. Ideally, we could
563 generate the new shake_vir, but test the veta value for convergence. This will take some thought. */
565 if (ir
->efep
!= efepNO
)
567 /* Sum up the foreign energy and dhdl terms for md and sd.
568 Currently done every step so that dhdl is correct in the .edr */
569 sum_dhdl(enerd
, state
->lambda
, ir
->fepvals
);
575 const bool bCalcEnerStep
= true;
576 upd_mdebin(mdebin
, doFreeEnergyPerturbation
, bCalcEnerStep
,
577 t
, mdatoms
->tmass
, enerd
, state
,
578 ir
->fepvals
, ir
->expandedvals
, state
->box
,
579 shake_vir
, force_vir
, total_vir
, pres
,
580 ekind
, mu_tot
, constr
);
582 const bool do_ene
= true;
583 const bool do_log
= true;
585 const bool do_dr
= ir
->nstdisreout
!= 0;
586 const bool do_or
= ir
->nstorireout
!= 0;
588 print_ebin(mdoutf_get_fp_ene(outf
), do_ene
, do_dr
, do_or
, do_log
? fplog
: nullptr,
590 eprNORMAL
, mdebin
, fcd
, groups
, &(ir
->opts
), awh
);
592 if (do_per_step(step
, ir
->nstlog
))
594 if (fflush(fplog
) != 0)
596 gmx_fatal(FARGS
, "Cannot flush logfile - maybe you are out of disk space?");
601 /* Print the remaining wall clock time for the run */
602 if (isMasterSimMasterRank(ms
, cr
) &&
603 (mdrunOptions
.verbose
|| gmx_got_usr_signal()))
607 fprintf(stderr
, "\n");
609 print_time(stderr
, walltime_accounting
, step
, ir
, cr
);
612 cycles
= wallcycle_stop(wcycle
, ewcSTEP
);
613 if (DOMAINDECOMP(cr
) && wcycle
)
615 dd_cycles_add(cr
->dd
, cycles
, ddCyclStep
);
618 /* increase the MD step number */
622 /* End of main MD loop */
624 /* Closing TNG files can include compressing data. Therefore it is good to do that
625 * before stopping the time measurements. */
626 mdoutf_tng_close(outf
);
628 /* Stop measuring walltime */
629 walltime_accounting_end_time(walltime_accounting
);
633 mimicCommunicator
.finalize();
636 if (!thisRankHasDuty(cr
, DUTY_PME
))
638 /* Tell the PME only node to finish */
639 gmx_pme_send_finish(cr
);
645 done_shellfc(fplog
, shellfc
, step_rel
);
647 // Clean up swapcoords
648 if (ir
->eSwapCoords
!= eswapNO
)
650 finish_swapcoords(ir
->swap
);
653 walltime_accounting_set_nsteps_done(walltime_accounting
, step_rel
);
655 destroy_enerdata(enerd
);