2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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/topology/atoms.h"
56 #include "energyelement.h"
57 #include "statepropagatordata.h"
61 struct gmx_multisim_t;
67 bool ShellFCElement::doShellsOrFlexConstraints(
68 const gmx_mtop_t *mtop, int nflexcon)
74 std::array<int, eptNR> n = countPtypes(nullptr, mtop);
75 return n[eptShell] != 0;
78 ShellFCElement::ShellFCElement(
79 StatePropagatorData *statePropagatorData,
80 EnergyElement *energyElement,
85 const t_inputrec *inputrec,
86 const MDAtoms *mdAtoms,
90 gmx_wallcycle *wcycle,
91 MdScheduleWorkload *mdScheduleWork,
93 ImdSession *imdSession,
96 const gmx_mtop_t *globalTopology) :
98 nextEnergyCalculationStep_(-1),
99 nextVirialCalculationStep_(-1),
100 nextFreeEnergyCalculationStep_(-1),
101 statePropagatorData_(statePropagatorData),
102 energyElement_(energyElement),
103 localTopology_(nullptr),
104 isDynamicBox_(isDynamicBox),
105 isVerbose_(isVerbose),
107 ddBalanceRegionHandler_(cr),
116 imdSession_(imdSession),
117 pull_work_(pull_work),
119 mdScheduleWork_(mdScheduleWork),
122 shellfc_ = init_shell_flexcon(
123 fplog, globalTopology,
124 constr_ ? constr_->numFlexibleConstraints() : 0,
125 inputrec->nstcalcenergy, DOMAINDECOMP(cr));
127 GMX_ASSERT(shellfc_, "ShellFCElement built, but init_shell_flexcon returned a nullptr");
129 if (!DOMAINDECOMP(cr))
131 // This was done in mdAlgorithmsSetupAtomData(), but shellfc
132 // won't be available outside this element.
133 make_local_shells(cr, mdAtoms->mdatoms(), shellfc_);
137 void ShellFCElement::scheduleTask(
138 Step step, Time time,
139 const RegisterRunFunctionPtr ®isterRunFunction)
141 unsigned int flags = (
142 GMX_FORCE_STATECHANGED |
143 GMX_FORCE_ALLFORCES |
144 (isDynamicBox_ ? GMX_FORCE_DYNAMICBOX : 0) |
145 (nextVirialCalculationStep_ == step ? GMX_FORCE_VIRIAL : 0) |
146 (nextEnergyCalculationStep_ == step ? GMX_FORCE_ENERGY : 0) |
147 (nextFreeEnergyCalculationStep_ == step ? GMX_FORCE_DHDL : 0));
149 (*registerRunFunction)(
150 std::make_unique<SimulatorRunFunction>(
151 [this, step, time, flags]()
152 {run(step, time, flags); }));
156 void ShellFCElement::elementSetup()
158 GMX_ASSERT(localTopology_, "Setup called before local topology was set.");
162 void ShellFCElement::run(Step step, Time time, unsigned int flags)
164 // Disabled functionality
165 gmx_multisim_t *ms = nullptr;
166 gmx_enfrot *enforcedRotation = nullptr;
167 t_graph *graph = nullptr;
169 auto x = statePropagatorData_->positionsView();
170 auto v = statePropagatorData_->velocitiesView();
171 auto forces = statePropagatorData_->forcesView();
172 auto box = statePropagatorData_->constBox();
173 auto lambda = ArrayRef<real>(); // disabled
174 history_t *hist = nullptr; // disabled
176 tensor force_vir = {{0}};
177 relax_shell_flexcon(fplog_, cr_, ms, isVerbose_,
178 enforcedRotation, step,
179 inputrec_, imdSession_, pull_work_, step == nextNSStep_,
180 static_cast<int>(flags), localTopology_,
181 constr_, energyElement_->enerdata(), fcd_,
182 statePropagatorData_->localNumAtoms(),
183 x, v, box, lambda, hist, forces, force_vir,
184 mdAtoms_->mdatoms(), nrnb_, wcycle_, graph,
185 shellfc_, fr_, mdScheduleWork_, time,
186 energyElement_->muTot(), vsite_,
187 ddBalanceRegionHandler_);
188 energyElement_->addToForceVirial(force_vir, step);
191 void ShellFCElement::elementTeardown()
193 done_shellfc(fplog_, shellfc_, nSteps_);
196 void ShellFCElement::setTopology(const gmx_localtop_t *top)
198 localTopology_ = top;
201 SignallerCallbackPtr ShellFCElement::registerNSCallback()
203 return std::make_unique<SignallerCallback>(
204 [this](Step step, Time gmx_unused time)
205 {this->nextNSStep_ = step; });
208 SignallerCallbackPtr ShellFCElement::
209 registerEnergyCallback(EnergySignallerEvent event)
211 if (event == EnergySignallerEvent::energyCalculationStep)
213 return std::make_unique<SignallerCallback>(
214 [this](Step step, Time)
215 {nextEnergyCalculationStep_ = step; });
217 if (event == EnergySignallerEvent::virialCalculationStep)
219 return std::make_unique<SignallerCallback>(
220 [this](Step step, Time){nextVirialCalculationStep_ = step; });
222 if (event == EnergySignallerEvent::freeEnergyCalculationStep)
224 return std::make_unique<SignallerCallback>(
225 [this](Step step, Time){nextFreeEnergyCalculationStep_ = step; });