2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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.
36 * \brief Defines the global reduction element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "computeglobalselement.h"
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/domdec/localtopologychecker.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdlib/constr.h"
52 #include "gromacs/mdlib/md_support.h"
53 #include "gromacs/mdlib/mdatoms.h"
54 #include "gromacs/mdlib/stat.h"
55 #include "gromacs/mdlib/update.h"
56 #include "gromacs/mdlib/vcm.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/mdtypes/group.h"
59 #include "gromacs/mdtypes/inputrec.h"
60 #include "gromacs/mdtypes/mdatom.h"
61 #include "gromacs/topology/topology.h"
63 #include "freeenergyperturbationdata.h"
64 #include "modularsimulator.h"
65 #include "simulatoralgorithm.h"
69 template<ComputeGlobalsAlgorithm algorithm>
70 ComputeGlobalsElement<algorithm>::ComputeGlobalsElement(StatePropagatorData* statePropagatorData,
71 EnergyData* energyData,
72 FreeEnergyPerturbationData* freeEnergyPerturbationData,
73 SimulationSignals* signals,
76 const MDLogger& mdlog,
78 const t_inputrec* inputrec,
79 const MDAtoms* mdAtoms,
81 gmx_wallcycle* wcycle,
83 const gmx_mtop_t& global_top,
84 Constraints* constr) :
85 energyReductionStep_(-1),
86 virialReductionStep_(-1),
87 vvSchedulingStep_(-1),
88 doStopCM_(inputrec->comm_mode != ComRemovalAlgorithm::No),
89 nstcomm_(inputrec->nstcomm),
90 nstglobalcomm_(nstglobalcomm),
91 lastStep_(inputrec->nsteps + inputrec->init_step),
92 initStep_(inputrec->init_step),
93 nullSignaller_(std::make_unique<SimulationSignaller>(nullptr, nullptr, nullptr, false, false)),
94 statePropagatorData_(statePropagatorData),
95 energyData_(energyData),
96 localTopology_(nullptr),
97 freeEnergyPerturbationData_(freeEnergyPerturbationData),
98 vcm_(global_top.groups, *inputrec),
104 top_global_(global_top),
111 reportComRemovalInfo(fplog, vcm_);
112 gstat_ = global_stat_init(inputrec_);
115 template<ComputeGlobalsAlgorithm algorithm>
116 ComputeGlobalsElement<algorithm>::~ComputeGlobalsElement()
118 global_stat_destroy(gstat_);
121 template<ComputeGlobalsAlgorithm algorithm>
122 void ComputeGlobalsElement<algorithm>::elementSetup()
124 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
126 if (doStopCM_ && !inputrec_->bContinuation)
128 // To minimize communication, compute_globals computes the COM velocity
129 // and the kinetic energy for the velocities without COM motion removed.
130 // Thus to get the kinetic energy without the COM contribution, we need
131 // to call compute_globals twice.
133 compute(-1, CGLO_GSTAT | CGLO_STOPCM, nullSignaller_.get(), false, true);
135 auto v = statePropagatorData_->velocitiesView();
136 // At initialization, do not pass x with acceleration-correction mode
137 // to avoid (incorrect) correction of the initial coordinates.
138 auto x = vcm_.mode == ComRemovalAlgorithm::LinearAccelerationCorrection
139 ? ArrayRefWithPadding<RVec>()
140 : statePropagatorData_->positionsView();
141 process_and_stopcm_grp(
142 fplog_, &vcm_, *mdAtoms_->mdatoms(), x.unpaddedArrayRef(), v.unpaddedArrayRef());
143 inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
146 unsigned int cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
147 | (energyData_->hasReadEkinFromCheckpoint() ? CGLO_READEKIN : 0));
149 if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
151 cglo_flags |= CGLO_PRESSURE | CGLO_CONSTRAINT;
154 compute(-1, cglo_flags, nullSignaller_.get(), false, true);
156 // Calculate the initial half step temperature, and save the ekinh_old
157 for (int i = 0; (i < inputrec_->opts.ngtc); i++)
159 copy_mat(energyData_->ekindata()->tcstat[i].ekinh, energyData_->ekindata()->tcstat[i].ekinh_old);
163 template<ComputeGlobalsAlgorithm algorithm>
164 void ComputeGlobalsElement<algorithm>::scheduleTask(Step step,
165 Time gmx_unused time,
166 const RegisterRunFunction& registerRunFunction)
168 const bool needComReduction = doStopCM_ && do_per_step(step, nstcomm_);
169 const bool needGlobalReduction = step == energyReductionStep_ || step == virialReductionStep_
170 || needComReduction || do_per_step(step, nstglobalcomm_);
172 // TODO: CGLO_GSTAT is only used for needToSumEkinhOld_, i.e. to signal that we do or do not
173 // sum the previous kinetic energy. We should simplify / clarify this.
175 if (algorithm == ComputeGlobalsAlgorithm::LeapFrog)
177 // With Leap-Frog we can skip compute_globals at
178 // non-communication steps, but we need to calculate
179 // the kinetic energy one step before communication.
181 // With leap-frog we also need to compute the half-step
182 // kinetic energy at the step before we need to write
183 // the full-step kinetic energy
184 const bool needEkinAtNextStep = (do_per_step(step + 1, nstglobalcomm_) || step + 1 == lastStep_);
186 if (!needGlobalReduction && !needEkinAtNextStep)
191 const bool doEnergy = step == energyReductionStep_;
192 int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
193 | (needComReduction ? CGLO_STOPCM : 0) | CGLO_TEMPERATURE | CGLO_PRESSURE
196 // Since we're already communicating at this step, we
197 // can propagate intra-simulation signals. Note that
198 // check_nstglobalcomm has the responsibility for
199 // choosing the value of nstglobalcomm which satisfies
200 // the need of the different signallers.
201 const bool doIntraSimSignal = true;
202 // Disable functionality
203 const bool doInterSimSignal = false;
205 // Make signaller to signal stop / reset / checkpointing signals
206 auto signaller = std::make_shared<SimulationSignaller>(
207 signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
209 registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
210 compute(step, flags, signaller.get(), true);
213 else if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
215 // For VV, we schedule two calls to compute globals per step.
216 if (step != vvSchedulingStep_)
218 // This is the first scheduling call for this step (positions & velocities at full time
219 // step) Set this as the current scheduling step
220 vvSchedulingStep_ = step;
222 // For vv, the state at the beginning of the step is positions at time t, velocities at time t - dt/2
223 // The first velocity propagation (+dt/2) therefore actually corresponds to the previous step.
224 // So we need information from the last step in the first half of the integration
225 if (!needGlobalReduction && !do_per_step(step - 1, nstglobalcomm_))
230 const bool doTemperature = step != initStep_ || inputrec_->bContinuation;
231 const bool doEnergy = step == energyReductionStep_;
233 int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
234 | (doTemperature ? CGLO_TEMPERATURE : 0) | CGLO_PRESSURE | CGLO_CONSTRAINT
235 | (needComReduction ? CGLO_STOPCM : 0) | CGLO_SCALEEKIN;
238 [this, step, flags]() { compute(step, flags, nullSignaller_.get(), false); });
242 // second call to compute_globals for this step
243 // Reset the scheduling step to avoid confusion if scheduling needs
244 // to be repeated (in case of unexpected simulation termination)
245 vvSchedulingStep_ = -1;
247 if (!needGlobalReduction)
251 int flags = CGLO_GSTAT | CGLO_CONSTRAINT;
253 // Since we're already communicating at this step, we
254 // can propagate intra-simulation signals. Note that
255 // check_nstglobalcomm has the responsibility for
256 // choosing the value of nstglobalcomm which satisfies
257 // the need of the different signallers.
258 const bool doIntraSimSignal = true;
259 // Disable functionality
260 const bool doInterSimSignal = false;
262 auto signaller = std::make_shared<SimulationSignaller>(
263 signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
265 registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
266 compute(step, flags, signaller.get(), true);
272 template<ComputeGlobalsAlgorithm algorithm>
273 void ComputeGlobalsElement<algorithm>::compute(gmx::Step step,
275 SimulationSignaller* signaller,
279 auto x = statePropagatorData_->positionsView().unpaddedArrayRef();
280 auto v = statePropagatorData_->velocitiesView().unpaddedArrayRef();
281 const auto* box = statePropagatorData_->constBox();
282 const auto* lastbox = useLastBox ? statePropagatorData_->constPreviousBox()
283 : statePropagatorData_->constBox();
285 if (DOMAINDECOMP(cr_) && dd_localTopologyChecker(*cr_->dd).shouldCheckNumberOfBondedInteractions())
287 flags |= CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS;
289 compute_globals(gstat_,
293 energyData_->ekindata(),
300 step != -1 ? wcycle_ : nullptr,
301 energyData_->enerdata(),
302 energyData_->forceVirial(step),
303 energyData_->constraintVirial(step),
304 energyData_->totalVirial(step),
305 energyData_->pressure(step),
306 (((flags & CGLO_ENERGY) != 0) && constr_ != nullptr) ? constr_->rmsdData()
307 : gmx::ArrayRef<real>{},
310 energyData_->needToSumEkinhOld(),
312 if (DOMAINDECOMP(cr_))
314 dd_localTopologyChecker(cr_->dd)->checkNumberOfBondedInteractions(localTopology_, x, box);
316 if (flags & CGLO_STOPCM && !isInit)
318 process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), x, v);
319 inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
323 template<ComputeGlobalsAlgorithm algorithm>
324 void ComputeGlobalsElement<algorithm>::setTopology(const gmx_localtop_t* top)
326 localTopology_ = top;
329 template<ComputeGlobalsAlgorithm algorithm>
330 std::optional<SignallerCallback> ComputeGlobalsElement<algorithm>::registerEnergyCallback(EnergySignallerEvent event)
332 if (event == EnergySignallerEvent::EnergyCalculationStep)
334 return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
336 if (event == EnergySignallerEvent::VirialCalculationStep)
338 return [this](Step step, Time /*unused*/) { virialReductionStep_ = step; };
343 template<ComputeGlobalsAlgorithm algorithm>
344 std::optional<SignallerCallback>
345 ComputeGlobalsElement<algorithm>::registerTrajectorySignallerCallback(TrajectoryEvent event)
347 if (event == TrajectoryEvent::EnergyWritingStep)
349 return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
354 //! Explicit template instantiation
356 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>;
357 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>;
361 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>::getElementPointerImpl(
362 LegacySimulatorData* legacySimulatorData,
363 ModularSimulatorAlgorithmBuilderHelper* builderHelper,
364 StatePropagatorData* statePropagatorData,
365 EnergyData* energyData,
366 FreeEnergyPerturbationData* freeEnergyPerturbationData,
367 GlobalCommunicationHelper* globalCommunicationHelper)
369 auto* element = builderHelper->storeElement(
370 std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>(
373 freeEnergyPerturbationData,
374 globalCommunicationHelper->simulationSignals(),
375 globalCommunicationHelper->nstglobalcomm(),
376 legacySimulatorData->fplog,
377 legacySimulatorData->mdlog,
378 legacySimulatorData->cr,
379 legacySimulatorData->inputrec,
380 legacySimulatorData->mdAtoms,
381 legacySimulatorData->nrnb,
382 legacySimulatorData->wcycle,
383 legacySimulatorData->fr,
384 legacySimulatorData->top_global,
385 legacySimulatorData->constr));
391 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>::getElementPointerImpl(
392 LegacySimulatorData* simulator,
393 ModularSimulatorAlgorithmBuilderHelper* builderHelper,
394 StatePropagatorData* statePropagatorData,
395 EnergyData* energyData,
396 FreeEnergyPerturbationData* freeEnergyPerturbationData,
397 GlobalCommunicationHelper* globalCommunicationHelper)
399 // We allow this element to be added multiple times to the call list, but we only want one
400 // actual element built
401 static const std::string key("vvComputeGlobalsElement");
403 const std::optional<std::any> cachedValue = builderHelper->builderData(key);
407 return std::any_cast<ISimulatorElement*>(cachedValue.value());
411 ISimulatorElement* vvComputeGlobalsElement = builderHelper->storeElement(
412 std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>(
415 freeEnergyPerturbationData,
416 globalCommunicationHelper->simulationSignals(),
417 globalCommunicationHelper->nstglobalcomm(),
426 simulator->top_global,
428 builderHelper->storeBuilderData(key, vvComputeGlobalsElement);
429 return vvComputeGlobalsElement;