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 shell / flex constraints element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "shellfcelement.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/topology/atoms.h"
56 #include "gromacs/topology/mtop_util.h"
58 #include "energyelement.h"
59 #include "freeenergyperturbationelement.h"
60 #include "statepropagatordata.h"
64 struct gmx_multisim_t;
70 bool ShellFCElement::doShellsOrFlexConstraints(const gmx_mtop_t& mtop, int nflexcon)
76 std::array<int, eptNR> n = gmx_mtop_particletype_count(mtop);
77 return n[eptShell] != 0;
80 ShellFCElement::ShellFCElement(StatePropagatorData* statePropagatorData,
81 EnergyElement* energyElement,
82 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
87 const t_inputrec* inputrec,
88 const MDAtoms* mdAtoms,
92 gmx_wallcycle* wcycle,
93 MdrunScheduleWorkload* runScheduleWork,
95 ImdSession* imdSession,
98 const gmx_mtop_t* globalTopology,
99 gmx_enfrot* enforcedRotation) :
101 nextEnergyCalculationStep_(-1),
102 nextVirialCalculationStep_(-1),
103 nextFreeEnergyCalculationStep_(-1),
104 statePropagatorData_(statePropagatorData),
105 energyElement_(energyElement),
106 freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
107 localTopology_(nullptr),
108 isDynamicBox_(isDynamicBox),
109 isVerbose_(isVerbose),
111 ddBalanceRegionHandler_(cr),
120 imdSession_(imdSession),
121 pull_work_(pull_work),
123 runScheduleWork_(runScheduleWork),
125 enforcedRotation_(enforcedRotation)
129 shellfc_ = init_shell_flexcon(fplog, globalTopology, constr_ ? constr_->numFlexibleConstraints() : 0,
130 inputrec->nstcalcenergy, DOMAINDECOMP(cr));
132 GMX_ASSERT(shellfc_, "ShellFCElement built, but init_shell_flexcon returned a nullptr");
134 if (!DOMAINDECOMP(cr))
136 // This was done in mdAlgorithmsSetupAtomData(), but shellfc
137 // won't be available outside this element.
138 make_local_shells(cr, mdAtoms->mdatoms(), shellfc_);
142 void ShellFCElement::scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction)
145 (GMX_FORCE_STATECHANGED | GMX_FORCE_ALLFORCES | (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0)
146 | (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0)
147 | (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0)
148 | (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0));
150 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
151 [this, step, time, flags]() { run(step, time, flags); }));
155 void ShellFCElement::elementSetup()
157 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
160 void ShellFCElement::run(Step step, Time time, unsigned int flags)
162 // Disabled functionality
163 gmx_multisim_t* ms = nullptr;
164 t_graph* graph = nullptr;
166 auto x = statePropagatorData_->positionsView();
167 auto v = statePropagatorData_->velocitiesView();
168 auto forces = statePropagatorData_->forcesView();
169 auto box = statePropagatorData_->constBox();
170 history_t* hist = nullptr; // disabled
172 tensor force_vir = { { 0 } };
173 // TODO: Make lambda const (needs some adjustments in lower force routines)
174 ArrayRef<real> lambda =
175 freeEnergyPerturbationElement_ ? freeEnergyPerturbationElement_->lambdaView() : lambda_;
177 fplog_, cr_, ms, isVerbose_, enforcedRotation_, step, inputrec_, imdSession_,
178 pull_work_, step == nextNSStep_, static_cast<int>(flags), localTopology_, constr_,
179 energyElement_->enerdata(), fcd_, statePropagatorData_->localNumAtoms(), x, v, box,
180 lambda, hist, forces, force_vir, mdAtoms_->mdatoms(), nrnb_, wcycle_, graph, shellfc_,
181 fr_, runScheduleWork_, time, energyElement_->muTot(), vsite_, ddBalanceRegionHandler_);
182 energyElement_->addToForceVirial(force_vir, step);
185 void ShellFCElement::elementTeardown()
187 done_shellfc(fplog_, shellfc_, nSteps_);
190 void ShellFCElement::setTopology(const gmx_localtop_t* top)
192 localTopology_ = top;
195 SignallerCallbackPtr ShellFCElement::registerNSCallback()
197 return std::make_unique<SignallerCallback>(
198 [this](Step step, Time gmx_unused time) { this->nextNSStep_ = step; });
201 SignallerCallbackPtr ShellFCElement::registerEnergyCallback(EnergySignallerEvent event)
203 if (event == EnergySignallerEvent::EnergyCalculationStep)
205 return std::make_unique<SignallerCallback>(
206 [this](Step step, Time /*unused*/) { nextEnergyCalculationStep_ = step; });
208 if (event == EnergySignallerEvent::VirialCalculationStep)
210 return std::make_unique<SignallerCallback>(
211 [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; });
213 if (event == EnergySignallerEvent::FreeEnergyCalculationStep)
215 return std::make_unique<SignallerCallback>(
216 [this](Step step, Time /*unused*/) { nextFreeEnergyCalculationStep_ = step; });