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 force element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "forceelement.h"
46 #include "gromacs/domdec/mdsetup.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/mdlib/constr.h"
49 #include "gromacs/mdlib/force.h"
50 #include "gromacs/mdlib/force_flags.h"
51 #include "gromacs/mdlib/mdatoms.h"
52 #include "gromacs/mdrun/shellfc.h"
53 #include "gromacs/mdtypes/forcebuffers.h"
54 #include "gromacs/mdtypes/forcerec.h"
55 #include "gromacs/mdtypes/inputrec.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/mdatom.h"
58 #include "gromacs/mdtypes/mdrunoptions.h"
59 #include "gromacs/mdtypes/simulation_workload.h"
60 #include "gromacs/pbcutil/pbc.h"
62 #include "energydata.h"
63 #include "freeenergyperturbationdata.h"
64 #include "modularsimulator.h"
65 #include "simulatoralgorithm.h"
66 #include "statepropagatordata.h"
70 struct gmx_multisim_t;
75 ForceElement::ForceElement(StatePropagatorData* statePropagatorData,
76 EnergyData* energyData,
77 FreeEnergyPerturbationData* freeEnergyPerturbationData,
82 const t_inputrec* inputrec,
83 const MDAtoms* mdAtoms,
86 gmx_wallcycle* wcycle,
87 MdrunScheduleWorkload* runScheduleWork,
88 VirtualSitesHandler* vsite,
89 ImdSession* imdSession,
92 const gmx_mtop_t& globalTopology,
93 gmx_enfrot* enforcedRotation) :
94 shellfc_(init_shell_flexcon(fplog,
96 constr ? constr->numFlexibleConstraints() : 0,
97 inputrec->nstcalcenergy,
98 haveDDAtomOrdering(*cr),
99 runScheduleWork->simulationWork.useGpuPme)),
100 doShellFC_(shellfc_ != nullptr),
102 nextEnergyCalculationStep_(-1),
103 nextVirialCalculationStep_(-1),
104 nextFreeEnergyCalculationStep_(-1),
105 statePropagatorData_(statePropagatorData),
106 energyData_(energyData),
107 freeEnergyPerturbationData_(freeEnergyPerturbationData),
108 localTopology_(nullptr),
109 isDynamicBox_(isDynamicBox),
110 isVerbose_(isVerbose),
111 nShellRelaxationSteps_(0),
112 ddBalanceRegionHandler_(cr),
113 longRangeNonbondeds_(std::make_unique<CpuPpLongRangeNonbondeds>(fr->n_tpi,
114 fr->ic->ewaldcoeff_q,
132 imdSession_(imdSession),
133 pull_work_(pull_work),
134 runScheduleWork_(runScheduleWork),
136 enforcedRotation_(enforcedRotation)
138 std::fill(lambda_.begin(), lambda_.end(), 0);
140 if (doShellFC_ && !haveDDAtomOrdering(*cr))
142 // This was done in mdAlgorithmsSetupAtomData(), but shellfc
143 // won't be available outside this element.
144 make_local_shells(cr, *mdAtoms->mdatoms(), shellfc_);
148 void ForceElement::scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction)
151 (GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0)
152 | (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0)
153 | (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0)
154 | (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0)
155 | (!doShellFC_ && nextNSStep_ == step ? GMX_FORCE_NS : 0));
157 registerRunFunction([this, step, time, flags]() {
160 run<true>(step, time, flags);
164 run<false>(step, time, flags);
169 void ForceElement::elementSetup()
171 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
174 template<bool doShellFC>
175 void ForceElement::run(Step step, Time time, unsigned int flags)
177 // Disabled functionality
178 gmx_multisim_t* ms = nullptr;
181 if (!haveDDAtomOrdering(*cr_) && (flags & GMX_FORCE_NS) && inputrecDynamicBox(inputrec_))
183 // TODO: Correcting the box is done in DomDecHelper (if using DD) or here (non-DD simulations).
184 // Think about unifying this responsibility, could this be done in one place?
185 auto* box = statePropagatorData_->box();
186 correct_box(fplog_, step, box);
189 /* The coordinates (x) are shifted (to get whole molecules)
191 * This is parallelized as well, and does communication too.
192 * Check comments in sim_util.c
194 auto x = statePropagatorData_->positionsView();
195 auto& forces = statePropagatorData_->forcesView();
196 const auto* box = statePropagatorData_->constBox();
197 history_t* hist = nullptr; // disabled
199 tensor force_vir = { { 0 } };
200 // TODO: Make lambda const (needs some adjustments in lower force routines)
201 ArrayRef<real> lambda =
202 freeEnergyPerturbationData_ ? freeEnergyPerturbationData_->lambdaView() : lambda_;
204 longRangeNonbondeds_->updateAfterPartition(*mdAtoms_->mdatoms());
208 auto v = statePropagatorData_->velocitiesView();
210 relax_shell_flexcon(fplog_,
220 static_cast<int>(flags),
223 energyData_->enerdata(),
224 statePropagatorData_->localNumAtoms(),
232 *mdAtoms_->mdatoms(),
233 longRangeNonbondeds_.get(),
240 energyData_->muTot(),
242 ddBalanceRegionHandler_);
243 nShellRelaxationSteps_++;
247 // Disabled functionality
249 gmx_edsam* ed = nullptr;
269 energyData_->enerdata(),
274 energyData_->muTot(),
277 longRangeNonbondeds_.get(),
278 static_cast<int>(flags),
279 ddBalanceRegionHandler_);
281 energyData_->addToForceVirial(force_vir, step);
284 void ForceElement::elementTeardown()
288 done_shellfc(fplog_, shellfc_, nShellRelaxationSteps_);
292 void ForceElement::setTopology(const gmx_localtop_t* top)
294 localTopology_ = top;
297 std::optional<SignallerCallback> ForceElement::registerNSCallback()
299 return [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; };
302 std::optional<SignallerCallback> ForceElement::registerEnergyCallback(EnergySignallerEvent event)
304 if (event == EnergySignallerEvent::EnergyCalculationStep)
306 return [this](Step step, Time /*unused*/) { nextEnergyCalculationStep_ = step; };
308 if (event == EnergySignallerEvent::VirialCalculationStep)
310 return [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; };
312 if (event == EnergySignallerEvent::FreeEnergyCalculationStep)
314 return [this](Step step, Time /*unused*/) { nextFreeEnergyCalculationStep_ = step; };
319 DomDecCallback ForceElement::registerDomDecCallback()
321 return [this]() { longRangeNonbondeds_->updateAfterPartition(*mdAtoms_->mdatoms()); };
325 ForceElement::getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
326 ModularSimulatorAlgorithmBuilderHelper* builderHelper,
327 StatePropagatorData* statePropagatorData,
328 EnergyData* energyData,
329 FreeEnergyPerturbationData* freeEnergyPerturbationData,
330 GlobalCommunicationHelper gmx_unused* globalCommunicationHelper,
331 ObservablesReducer* /*observablesReducer*/)
333 const bool isVerbose = legacySimulatorData->mdrunOptions.verbose;
334 const bool isDynamicBox = inputrecDynamicBox(legacySimulatorData->inputrec);
335 return builderHelper->storeElement(
336 std::make_unique<ForceElement>(statePropagatorData,
338 freeEnergyPerturbationData,
341 legacySimulatorData->fplog,
342 legacySimulatorData->cr,
343 legacySimulatorData->inputrec,
344 legacySimulatorData->mdAtoms,
345 legacySimulatorData->nrnb,
346 legacySimulatorData->fr,
347 legacySimulatorData->wcycle,
348 legacySimulatorData->runScheduleWork,
349 legacySimulatorData->vsite,
350 legacySimulatorData->imdSession,
351 legacySimulatorData->pull_work,
352 legacySimulatorData->constr,
353 legacySimulatorData->top_global,
354 legacySimulatorData->enforcedRotation));