2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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/inputrec.h"
54 #include "gromacs/mdtypes/mdatom.h"
55 #include "gromacs/pbcutil/pbc.h"
57 #include "energyelement.h"
58 #include "freeenergyperturbationelement.h"
59 #include "statepropagatordata.h"
63 struct gmx_multisim_t;
68 ForceElement::ForceElement(StatePropagatorData* statePropagatorData,
69 EnergyElement* energyElement,
70 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
75 const t_inputrec* inputrec,
76 const MDAtoms* mdAtoms,
80 gmx_wallcycle* wcycle,
81 MdrunScheduleWorkload* runScheduleWork,
82 VirtualSitesHandler* vsite,
83 ImdSession* imdSession,
86 const gmx_mtop_t* globalTopology,
87 gmx_enfrot* enforcedRotation) :
88 shellfc_(init_shell_flexcon(fplog,
90 constr ? constr->numFlexibleConstraints() : 0,
91 inputrec->nstcalcenergy,
93 doShellFC_(shellfc_ != nullptr),
95 nextEnergyCalculationStep_(-1),
96 nextVirialCalculationStep_(-1),
97 nextFreeEnergyCalculationStep_(-1),
98 statePropagatorData_(statePropagatorData),
99 energyElement_(energyElement),
100 freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
101 localTopology_(nullptr),
102 isDynamicBox_(isDynamicBox),
103 isVerbose_(isVerbose),
104 nShellRelaxationSteps_(0),
105 ddBalanceRegionHandler_(cr),
115 imdSession_(imdSession),
116 pull_work_(pull_work),
118 runScheduleWork_(runScheduleWork),
120 enforcedRotation_(enforcedRotation)
124 if (doShellFC_ && !DOMAINDECOMP(cr))
126 // This was done in mdAlgorithmsSetupAtomData(), but shellfc
127 // won't be available outside this element.
128 make_local_shells(cr, mdAtoms->mdatoms(), shellfc_);
132 void ForceElement::scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction)
135 (GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0)
136 | (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0)
137 | (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0)
138 | (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0)
139 | (!doShellFC_ && nextNSStep_ == step ? GMX_FORCE_NS : 0));
141 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>([this, step, time, flags]() {
144 run<true>(step, time, flags);
148 run<false>(step, time, flags);
153 void ForceElement::elementSetup()
155 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
158 template<bool doShellFC>
159 void ForceElement::run(Step step, Time time, unsigned int flags)
161 // Disabled functionality
162 gmx_multisim_t* ms = nullptr;
165 if (!DOMAINDECOMP(cr_) && (flags & GMX_FORCE_NS) && inputrecDynamicBox(inputrec_))
167 // TODO: Correcting the box is done in DomDecHelper (if using DD) or here (non-DD simulations).
168 // Think about unifying this responsibility, could this be done in one place?
169 auto box = statePropagatorData_->box();
170 correct_box(fplog_, step, box);
173 /* The coordinates (x) are shifted (to get whole molecules)
175 * This is parallelized as well, and does communication too.
176 * Check comments in sim_util.c
178 auto x = statePropagatorData_->positionsView();
179 auto forces = statePropagatorData_->forcesView();
180 auto box = statePropagatorData_->constBox();
181 history_t* hist = nullptr; // disabled
183 tensor force_vir = { { 0 } };
184 // TODO: Make lambda const (needs some adjustments in lower force routines)
185 ArrayRef<real> lambda =
186 freeEnergyPerturbationElement_ ? freeEnergyPerturbationElement_->lambdaView() : lambda_;
190 auto v = statePropagatorData_->velocitiesView();
193 fplog_, cr_, ms, isVerbose_, enforcedRotation_, step, inputrec_, imdSession_,
194 pull_work_, step == nextNSStep_, static_cast<int>(flags), localTopology_, constr_,
195 energyElement_->enerdata(), fcd_, statePropagatorData_->localNumAtoms(), x, v, box,
196 lambda, hist, forces, force_vir, mdAtoms_->mdatoms(), nrnb_, wcycle_, shellfc_, fr_,
197 runScheduleWork_, time, energyElement_->muTot(), vsite_, ddBalanceRegionHandler_);
198 nShellRelaxationSteps_++;
202 // Disabled functionality
204 gmx_edsam* ed = nullptr;
206 do_force(fplog_, cr_, ms, inputrec_, awh, enforcedRotation_, imdSession_, pull_work_, step,
207 nrnb_, wcycle_, localTopology_, box, x, hist, forces, force_vir, mdAtoms_->mdatoms(),
208 energyElement_->enerdata(), fcd_, lambda, fr_, runScheduleWork_, vsite_,
209 energyElement_->muTot(), time, ed, static_cast<int>(flags), ddBalanceRegionHandler_);
211 energyElement_->addToForceVirial(force_vir, step);
214 void ForceElement::elementTeardown()
218 done_shellfc(fplog_, shellfc_, nShellRelaxationSteps_);
222 void ForceElement::setTopology(const gmx_localtop_t* top)
224 localTopology_ = top;
227 SignallerCallbackPtr ForceElement::registerNSCallback()
229 return std::make_unique<SignallerCallback>(
230 [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; });
233 SignallerCallbackPtr ForceElement::registerEnergyCallback(EnergySignallerEvent event)
235 if (event == EnergySignallerEvent::EnergyCalculationStep)
237 return std::make_unique<SignallerCallback>(
238 [this](Step step, Time /*unused*/) { nextEnergyCalculationStep_ = step; });
240 if (event == EnergySignallerEvent::VirialCalculationStep)
242 return std::make_unique<SignallerCallback>(
243 [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; });
245 if (event == EnergySignallerEvent::FreeEnergyCalculationStep)
247 return std::make_unique<SignallerCallback>(
248 [this](Step step, Time /*unused*/) { nextFreeEnergyCalculationStep_ = step; });