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/partition.h"
47 #include "gromacs/gmxlib/network.h"
48 #include "gromacs/gmxlib/nrnb.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/constr.h"
51 #include "gromacs/mdlib/md_support.h"
52 #include "gromacs/mdlib/mdatoms.h"
53 #include "gromacs/mdlib/stat.h"
54 #include "gromacs/mdlib/update.h"
55 #include "gromacs/mdtypes/commrec.h"
56 #include "gromacs/mdtypes/group.h"
57 #include "gromacs/mdtypes/inputrec.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/mdtypes/mdatom.h"
60 #include "gromacs/topology/topology.h"
62 #include "freeenergyperturbationdata.h"
63 #include "modularsimulator.h"
64 #include "simulatoralgorithm.h"
68 template<ComputeGlobalsAlgorithm algorithm>
69 ComputeGlobalsElement<algorithm>::ComputeGlobalsElement(StatePropagatorData* statePropagatorData,
70 EnergyData* energyData,
71 FreeEnergyPerturbationData* freeEnergyPerturbationData,
72 SimulationSignals* signals,
75 const MDLogger& mdlog,
77 const t_inputrec* inputrec,
78 const MDAtoms* mdAtoms,
80 gmx_wallcycle* wcycle,
82 const gmx_mtop_t& global_top,
83 Constraints* constr) :
84 energyReductionStep_(-1),
85 virialReductionStep_(-1),
86 vvSchedulingStep_(-1),
87 doStopCM_(inputrec->comm_mode != ComRemovalAlgorithm::No),
88 nstcomm_(inputrec->nstcomm),
89 nstglobalcomm_(nstglobalcomm),
90 lastStep_(inputrec->nsteps + inputrec->init_step),
91 initStep_(inputrec->init_step),
92 nullSignaller_(std::make_unique<SimulationSignaller>(nullptr, nullptr, nullptr, false, false)),
93 totalNumberOfBondedInteractions_(0),
94 shouldCheckNumberOfBondedInteractions_(false),
95 statePropagatorData_(statePropagatorData),
96 energyData_(energyData),
97 localTopology_(nullptr),
98 freeEnergyPerturbationData_(freeEnergyPerturbationData),
99 vcm_(global_top.groups, *inputrec),
105 top_global_(global_top),
112 reportComRemovalInfo(fplog, vcm_);
113 gstat_ = global_stat_init(inputrec_);
116 template<ComputeGlobalsAlgorithm algorithm>
117 ComputeGlobalsElement<algorithm>::~ComputeGlobalsElement()
119 global_stat_destroy(gstat_);
122 template<ComputeGlobalsAlgorithm algorithm>
123 void ComputeGlobalsElement<algorithm>::elementSetup()
125 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
127 if (doStopCM_ && !inputrec_->bContinuation)
129 // To minimize communication, compute_globals computes the COM velocity
130 // and the kinetic energy for the velocities without COM motion removed.
131 // Thus to get the kinetic energy without the COM contribution, we need
132 // to call compute_globals twice.
134 compute(-1, CGLO_GSTAT | CGLO_STOPCM, nullSignaller_.get(), false, true);
136 auto v = statePropagatorData_->velocitiesView();
137 // At initialization, do not pass x with acceleration-correction mode
138 // to avoid (incorrect) correction of the initial coordinates.
139 auto x = vcm_.mode == ComRemovalAlgorithm::LinearAccelerationCorrection
140 ? ArrayRefWithPadding<RVec>()
141 : statePropagatorData_->positionsView();
142 process_and_stopcm_grp(
143 fplog_, &vcm_, *mdAtoms_->mdatoms(), x.unpaddedArrayRef(), v.unpaddedArrayRef());
144 inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
147 unsigned int cglo_flags = (CGLO_TEMPERATURE | CGLO_GSTAT
148 | (energyData_->hasReadEkinFromCheckpoint() ? CGLO_READEKIN : 0));
150 if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
152 cglo_flags |= CGLO_PRESSURE | CGLO_CONSTRAINT;
155 compute(-1, cglo_flags, nullSignaller_.get(), false, true);
157 // Calculate the initial half step temperature, and save the ekinh_old
158 for (int i = 0; (i < inputrec_->opts.ngtc); i++)
160 copy_mat(energyData_->ekindata()->tcstat[i].ekinh, energyData_->ekindata()->tcstat[i].ekinh_old);
164 template<ComputeGlobalsAlgorithm algorithm>
165 void ComputeGlobalsElement<algorithm>::scheduleTask(Step step,
166 Time gmx_unused time,
167 const RegisterRunFunction& registerRunFunction)
169 const bool needComReduction = doStopCM_ && do_per_step(step, nstcomm_);
170 const bool needGlobalReduction = step == energyReductionStep_ || step == virialReductionStep_
171 || needComReduction || do_per_step(step, nstglobalcomm_);
173 // TODO: CGLO_GSTAT is only used for needToSumEkinhOld_, i.e. to signal that we do or do not
174 // sum the previous kinetic energy. We should simplify / clarify this.
176 if (algorithm == ComputeGlobalsAlgorithm::LeapFrog)
178 // With Leap-Frog we can skip compute_globals at
179 // non-communication steps, but we need to calculate
180 // the kinetic energy one step before communication.
182 // With leap-frog we also need to compute the half-step
183 // kinetic energy at the step before we need to write
184 // the full-step kinetic energy
185 const bool needEkinAtNextStep = (do_per_step(step + 1, nstglobalcomm_) || step + 1 == lastStep_);
187 if (!needGlobalReduction && !needEkinAtNextStep)
192 const bool doEnergy = step == energyReductionStep_;
193 int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
194 | (needComReduction ? CGLO_STOPCM : 0) | CGLO_TEMPERATURE | CGLO_PRESSURE
197 // Since we're already communicating at this step, we
198 // can propagate intra-simulation signals. Note that
199 // check_nstglobalcomm has the responsibility for
200 // choosing the value of nstglobalcomm which satisfies
201 // the need of the different signallers.
202 const bool doIntraSimSignal = true;
203 // Disable functionality
204 const bool doInterSimSignal = false;
206 // Make signaller to signal stop / reset / checkpointing signals
207 auto signaller = std::make_shared<SimulationSignaller>(
208 signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
210 registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
211 compute(step, flags, signaller.get(), true);
214 else if (algorithm == ComputeGlobalsAlgorithm::VelocityVerlet)
216 // For VV, we schedule two calls to compute globals per step.
217 if (step != vvSchedulingStep_)
219 // This is the first scheduling call for this step (positions & velocities at full time
220 // step) Set this as the current scheduling step
221 vvSchedulingStep_ = step;
223 // For vv, the state at the beginning of the step is positions at time t, velocities at time t - dt/2
224 // The first velocity propagation (+dt/2) therefore actually corresponds to the previous step.
225 // So we need information from the last step in the first half of the integration
226 if (!needGlobalReduction && !do_per_step(step - 1, nstglobalcomm_))
231 const bool doTemperature = step != initStep_ || inputrec_->bContinuation;
232 const bool doEnergy = step == energyReductionStep_;
234 int flags = (needGlobalReduction ? CGLO_GSTAT : 0) | (doEnergy ? CGLO_ENERGY : 0)
235 | (doTemperature ? CGLO_TEMPERATURE : 0) | CGLO_PRESSURE | CGLO_CONSTRAINT
236 | (needComReduction ? CGLO_STOPCM : 0) | CGLO_SCALEEKIN;
239 [this, step, flags]() { compute(step, flags, nullSignaller_.get(), false); });
243 // second call to compute_globals for this step
244 // Reset the scheduling step to avoid confusion if scheduling needs
245 // to be repeated (in case of unexpected simulation termination)
246 vvSchedulingStep_ = -1;
248 if (!needGlobalReduction)
252 int flags = CGLO_GSTAT | CGLO_CONSTRAINT;
254 // Since we're already communicating at this step, we
255 // can propagate intra-simulation signals. Note that
256 // check_nstglobalcomm has the responsibility for
257 // choosing the value of nstglobalcomm which satisfies
258 // the need of the different signallers.
259 const bool doIntraSimSignal = true;
260 // Disable functionality
261 const bool doInterSimSignal = false;
263 auto signaller = std::make_shared<SimulationSignaller>(
264 signals_, cr_, nullptr, doInterSimSignal, doIntraSimSignal);
266 registerRunFunction([this, step, flags, signaller = std::move(signaller)]() {
267 compute(step, flags, signaller.get(), true);
273 template<ComputeGlobalsAlgorithm algorithm>
274 void ComputeGlobalsElement<algorithm>::compute(gmx::Step step,
276 SimulationSignaller* signaller,
280 auto x = statePropagatorData_->positionsView().unpaddedArrayRef();
281 auto v = statePropagatorData_->velocitiesView().unpaddedArrayRef();
282 auto box = statePropagatorData_->constBox();
283 auto lastbox = useLastBox ? statePropagatorData_->constPreviousBox()
284 : statePropagatorData_->constBox();
291 energyData_->ekindata(),
298 step != -1 ? wcycle_ : nullptr,
299 energyData_->enerdata(),
300 energyData_->forceVirial(step),
301 energyData_->constraintVirial(step),
302 energyData_->totalVirial(step),
303 energyData_->pressure(step),
304 (((flags & CGLO_ENERGY) != 0) && constr_ != nullptr) ? constr_->rmsdData()
305 : gmx::ArrayRef<real>{},
308 &totalNumberOfBondedInteractions_,
309 energyData_->needToSumEkinhOld(),
310 flags | (shouldCheckNumberOfBondedInteractions_ ? CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS : 0));
311 checkNumberOfBondedInteractions(
312 mdlog_, cr_, totalNumberOfBondedInteractions_, top_global_, localTopology_, x, box, &shouldCheckNumberOfBondedInteractions_);
313 if (flags & CGLO_STOPCM && !isInit)
315 process_and_stopcm_grp(fplog_, &vcm_, *mdAtoms_->mdatoms(), x, v);
316 inc_nrnb(nrnb_, eNR_STOPCM, mdAtoms_->mdatoms()->homenr);
320 template<ComputeGlobalsAlgorithm algorithm>
321 CheckBondedInteractionsCallback ComputeGlobalsElement<algorithm>::getCheckNumberOfBondedInteractionsCallback()
323 return [this]() { needToCheckNumberOfBondedInteractions(); };
326 template<ComputeGlobalsAlgorithm algorithm>
327 void ComputeGlobalsElement<algorithm>::needToCheckNumberOfBondedInteractions()
329 shouldCheckNumberOfBondedInteractions_ = true;
332 template<ComputeGlobalsAlgorithm algorithm>
333 void ComputeGlobalsElement<algorithm>::setTopology(const gmx_localtop_t* top)
335 localTopology_ = top;
338 template<ComputeGlobalsAlgorithm algorithm>
339 std::optional<SignallerCallback> ComputeGlobalsElement<algorithm>::registerEnergyCallback(EnergySignallerEvent event)
341 if (event == EnergySignallerEvent::EnergyCalculationStep)
343 return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
345 if (event == EnergySignallerEvent::VirialCalculationStep)
347 return [this](Step step, Time /*unused*/) { virialReductionStep_ = step; };
352 template<ComputeGlobalsAlgorithm algorithm>
353 std::optional<SignallerCallback>
354 ComputeGlobalsElement<algorithm>::registerTrajectorySignallerCallback(TrajectoryEvent event)
356 if (event == TrajectoryEvent::EnergyWritingStep)
358 return [this](Step step, Time /*unused*/) { energyReductionStep_ = step; };
363 //! Explicit template instantiation
365 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>;
366 template class ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>;
370 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>::getElementPointerImpl(
371 LegacySimulatorData* legacySimulatorData,
372 ModularSimulatorAlgorithmBuilderHelper* builderHelper,
373 StatePropagatorData* statePropagatorData,
374 EnergyData* energyData,
375 FreeEnergyPerturbationData* freeEnergyPerturbationData,
376 GlobalCommunicationHelper* globalCommunicationHelper)
378 auto* element = builderHelper->storeElement(
379 std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>>(
382 freeEnergyPerturbationData,
383 globalCommunicationHelper->simulationSignals(),
384 globalCommunicationHelper->nstglobalcomm(),
385 legacySimulatorData->fplog,
386 legacySimulatorData->mdlog,
387 legacySimulatorData->cr,
388 legacySimulatorData->inputrec,
389 legacySimulatorData->mdAtoms,
390 legacySimulatorData->nrnb,
391 legacySimulatorData->wcycle,
392 legacySimulatorData->fr,
393 legacySimulatorData->top_global,
394 legacySimulatorData->constr));
396 // TODO: Remove this when DD can reduce bonded interactions independently (#3421)
397 auto* castedElement = static_cast<ComputeGlobalsElement<ComputeGlobalsAlgorithm::LeapFrog>*>(element);
398 globalCommunicationHelper->setCheckBondedInteractionsCallback(
399 castedElement->getCheckNumberOfBondedInteractionsCallback());
405 ISimulatorElement* ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>::getElementPointerImpl(
406 LegacySimulatorData* simulator,
407 ModularSimulatorAlgorithmBuilderHelper* builderHelper,
408 StatePropagatorData* statePropagatorData,
409 EnergyData* energyData,
410 FreeEnergyPerturbationData* freeEnergyPerturbationData,
411 GlobalCommunicationHelper* globalCommunicationHelper)
413 // We allow this element to be added multiple times to the call list, but we only want one
414 // actual element built
415 static const std::string key("vvComputeGlobalsElement");
417 const std::optional<std::any> cachedValue = builderHelper->getStoredValue(key);
421 return std::any_cast<ISimulatorElement*>(cachedValue.value());
425 ISimulatorElement* vvComputeGlobalsElement = builderHelper->storeElement(
426 std::make_unique<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>>(
429 freeEnergyPerturbationData,
430 globalCommunicationHelper->simulationSignals(),
431 globalCommunicationHelper->nstglobalcomm(),
440 simulator->top_global,
443 // TODO: Remove this when DD can reduce bonded interactions independently (#3421)
444 auto* castedElement =
445 static_cast<ComputeGlobalsElement<ComputeGlobalsAlgorithm::VelocityVerlet>*>(
446 vvComputeGlobalsElement);
447 globalCommunicationHelper->setCheckBondedInteractionsCallback(
448 castedElement->getCheckNumberOfBondedInteractionsCallback());
449 builderHelper->storeValue(key, vvComputeGlobalsElement);
450 return vvComputeGlobalsElement;