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 constraint element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
44 #include "constraintelement.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/mdtypes/state.h"
49 #include "gromacs/utility/fatalerror.h"
51 #include "energyelement.h"
52 #include "freeenergyperturbationelement.h"
53 #include "statepropagatordata.h"
57 template<ConstraintVariable variable>
58 ConstraintsElement<variable>::ConstraintsElement(Constraints* constr,
59 StatePropagatorData* statePropagatorData,
60 EnergyElement* energyElement,
61 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
64 const t_inputrec* inputrec,
65 const t_mdatoms* mdAtoms) :
66 nextVirialCalculationStep_(-1),
67 nextEnergyWritingStep_(-1),
68 nextLogWritingStep_(-1),
69 isMasterRank_(isMaster),
70 statePropagatorData_(statePropagatorData),
71 energyElement_(energyElement),
72 freeEnergyPerturbationElement_(freeEnergyPerturbationElement),
78 GMX_ASSERT(constr_, "Constraint element created but constr == nullptr");
81 template<ConstraintVariable variable>
82 void ConstraintsElement<variable>::elementSetup()
84 if (!inputrec_->bContinuation
85 && ((variable == ConstraintVariable::Positions && inputrec_->eI == eiMD)
86 || (variable == ConstraintVariable::Velocities && inputrec_->eI == eiVV)))
88 const real lambdaBonded = freeEnergyPerturbationElement_
89 ? freeEnergyPerturbationElement_->constLambdaView()[efptBONDED]
91 // Constrain the initial coordinates and velocities
92 do_constrain_first(fplog_, constr_, inputrec_, mdAtoms_, statePropagatorData_->localNumAtoms(),
93 statePropagatorData_->positionsView(), statePropagatorData_->velocitiesView(),
94 statePropagatorData_->box(), lambdaBonded);
98 if (inputrec_->eConstrAlg == econtLINCS)
100 fprintf(fplog_, "RMS relative constraint deviation after constraining: %.2e\n",
107 template<ConstraintVariable variable>
108 void ConstraintsElement<variable>::scheduleTask(Step step,
109 Time gmx_unused time,
110 const RegisterRunFunctionPtr& registerRunFunction)
112 bool calculateVirial = (step == nextVirialCalculationStep_);
113 bool writeLog = (step == nextLogWritingStep_);
114 bool writeEnergy = (step == nextEnergyWritingStep_);
116 // register constraining
117 (*registerRunFunction)(std::make_unique<SimulatorRunFunction>(
118 [this, step, calculateVirial, writeLog, writeEnergy]() {
119 apply(step, calculateVirial, writeLog, writeEnergy);
123 template<ConstraintVariable variable>
124 void ConstraintsElement<variable>::apply(Step step, bool calculateVirial, bool writeLog, bool writeEnergy)
128 rvec *x, *xprime, *min_proj, *v;
130 // disabled functionality
132 real* dvdlambda = nullptr;
136 case ConstraintVariable::Positions:
137 x = as_rvec_array(statePropagatorData_->previousPositionsView().paddedArrayRef().data());
138 xprime = as_rvec_array(statePropagatorData_->positionsView().paddedArrayRef().data());
140 v = as_rvec_array(statePropagatorData_->velocitiesView().paddedArrayRef().data());
142 case ConstraintVariable::Velocities:
143 x = as_rvec_array(statePropagatorData_->positionsView().paddedArrayRef().data());
144 xprime = as_rvec_array(statePropagatorData_->velocitiesView().paddedArrayRef().data());
145 min_proj = as_rvec_array(statePropagatorData_->velocitiesView().paddedArrayRef().data());
148 default: gmx_fatal(FARGS, "Constraint algorithm not implemented for modular simulator.");
151 constr_->apply(writeLog, writeEnergy, step, 1, 1.0, x, xprime, min_proj, statePropagatorData_->box(),
152 lambda, dvdlambda, v, calculateVirial ? &vir_con : nullptr, variable);
156 if (inputrec_->eI == eiVV)
158 // For some reason, the shake virial in VV is reset twice a step.
159 // Energy element will only do this once per step.
160 // TODO: Investigate this
161 clear_mat(energyElement_->constraintVirial(step));
163 energyElement_->addToConstraintVirial(vir_con, step);
167 template<ConstraintVariable variable>
168 SignallerCallbackPtr ConstraintsElement<variable>::registerEnergyCallback(EnergySignallerEvent event)
170 if (event == EnergySignallerEvent::VirialCalculationStep)
172 return std::make_unique<SignallerCallback>(
173 [this](Step step, Time /*unused*/) { nextVirialCalculationStep_ = step; });
178 template<ConstraintVariable variable>
179 SignallerCallbackPtr ConstraintsElement<variable>::registerTrajectorySignallerCallback(TrajectoryEvent event)
181 if (event == TrajectoryEvent::EnergyWritingStep)
183 return std::make_unique<SignallerCallback>(
184 [this](Step step, Time /*unused*/) { nextEnergyWritingStep_ = step; });
189 template<ConstraintVariable variable>
190 SignallerCallbackPtr ConstraintsElement<variable>::registerLoggingCallback()
192 return std::make_unique<SignallerCallback>(
193 [this](Step step, Time /*unused*/) { nextLogWritingStep_ = step; });
196 //! Explicit template initialization
198 template class ConstraintsElement<ConstraintVariable::Positions>;
199 template class ConstraintsElement<ConstraintVariable::Velocities>;