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/math/vectypes.h"
47 #include "gromacs/mdlib/force.h"
48 #include "gromacs/mdlib/force_flags.h"
49 #include "gromacs/mdlib/mdatoms.h"
50 #include "gromacs/mdtypes/mdatom.h"
52 #include "energyelement.h"
53 #include "freeenergyperturbationelement.h"
54 #include "statepropagatordata.h"
58 struct gmx_multisim_t;
64 ForceElement::ForceElement(StatePropagatorData* statePropagatorData,
65 EnergyElement* energyElement,
66 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
70 const t_inputrec* inputrec,
71 const MDAtoms* mdAtoms,
75 gmx_wallcycle* wcycle,
76 MdrunScheduleWorkload* runScheduleWork,
78 ImdSession* imdSession,
80 gmx_enfrot* enforcedRotation) :
82 nextEnergyCalculationStep_(-1),
83 nextVirialCalculationStep_(-1),
84 nextFreeEnergyCalculationStep_(-1),
85 statePropagatorData_(statePropagatorData),
86 energyElement_(energyElement),
87 freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
88 localTopology_(nullptr),
89 isDynamicBox_(isDynamicBox),
90 ddBalanceRegionHandler_(cr),
100 imdSession_(imdSession),
101 pull_work_(pull_work),
103 runScheduleWork_(runScheduleWork),
104 enforcedRotation_(enforcedRotation)
109 void ForceElement::scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction)
112 (GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0)
113 | (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0)
114 | (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0)
115 | (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0)
116 | (nextNSStep_ == step ? GMX_FORCE_NS : 0));
118 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
119 [this, step, time, flags]() { run(step, time, flags); }));
122 void ForceElement::elementSetup()
124 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
127 void ForceElement::run(Step step, Time time, unsigned int flags)
129 // Disabled functionality
131 gmx_edsam* ed = nullptr;
132 gmx_multisim_t* ms = nullptr;
133 t_graph* graph = nullptr;
135 /* The coordinates (x) are shifted (to get whole molecules)
137 * This is parallelized as well, and does communication too.
138 * Check comments in sim_util.c
140 auto x = statePropagatorData_->positionsView();
141 auto forces = statePropagatorData_->forcesView();
142 auto box = statePropagatorData_->constBox();
143 history_t* hist = nullptr; // disabled
145 tensor force_vir = { { 0 } };
146 // TODO: Make lambda const (needs some adjustments in lower force routines)
147 ArrayRef<real> lambda =
148 freeEnergyPerturbationElement_ ? freeEnergyPerturbationElement_->lambdaView() : lambda_;
150 do_force(fplog_, cr_, ms, inputrec_, awh, enforcedRotation_, imdSession_, pull_work_, step,
151 nrnb_, wcycle_, localTopology_, box, x, hist, forces, force_vir, mdAtoms_->mdatoms(),
152 energyElement_->enerdata(), fcd_, lambda, graph, fr_, runScheduleWork_, vsite_,
153 energyElement_->muTot(), time, ed, static_cast<int>(flags), ddBalanceRegionHandler_);
154 energyElement_->addToForceVirial(force_vir, step);
157 void ForceElement::setTopology(const gmx_localtop_t* top)
159 localTopology_ = top;
162 SignallerCallbackPtr ForceElement::registerNSCallback()
164 return std::make_unique<SignallerCallback>(
165 [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; });
168 SignallerCallbackPtr ForceElement::registerEnergyCallback(EnergySignallerEvent event)
170 if (event == EnergySignallerEvent::EnergyCalculationStep)
172 return std::make_unique<SignallerCallback>(
173 [this](Step step, Time /*unused*/) { nextEnergyCalculationStep_ = step; });
175 if (event == EnergySignallerEvent::VirialCalculationStep)
177 return std::make_unique<SignallerCallback>(
178 [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; });
180 if (event == EnergySignallerEvent::FreeEnergyCalculationStep)
182 return std::make_unique<SignallerCallback>(
183 [this](Step step, Time /*unused*/) { nextFreeEnergyCalculationStep_ = step; });