Annotate modular simulator headers with exposure level
[gromacs.git] / src / gromacs / modularsimulator / propagator.h
blob3528fbd322925910d8ebcd8ca55d4f3da52b04ae
1 /*
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.
35 /*! \internal \file
36 * \brief Declares the propagator element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
41 * This header is only used within the modular simulator module
44 #ifndef GMX_MODULARSIMULATOR_PROPAGATOR_H
45 #define GMX_MODULARSIMULATOR_PROPAGATOR_H
47 #include <vector>
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/real.h"
53 #include "modularsimulatorinterfaces.h"
55 struct gmx_wallcycle;
57 namespace gmx
59 class MDAtoms;
60 class StatePropagatorData;
62 //! \addtogroup module_modularsimulator
63 //! \{
65 /*! \brief The different integration types we know about
67 * PositionsOnly:
68 * Moves the position vector by the given time step
69 * VelocitiesOnly:
70 * Moves the velocity vector by the given time step
71 * LeapFrog:
72 * Is a manual fusion of the previous two propagators
73 * VelocityVerletPositionsAndVelocities:
74 * Is a manual fusion of VelocitiesOnly and PositionsOnly,
75 * where VelocitiesOnly is only propagated by half the
76 * time step of the positions.
78 enum class IntegrationStep
80 PositionsOnly,
81 VelocitiesOnly,
82 LeapFrog,
83 VelocityVerletPositionsAndVelocities,
84 Count
87 //! Sets the number of different velocity scaling values
88 enum class NumVelocityScalingValues
90 None, //!< No velocity scaling (either this step or ever)
91 Single, //!< Single T-scaling value (either one group or all values =1)
92 Multiple, //!< Multiple T-scaling values, need to use T-group indices
93 Count
96 //! Sets the type of Parrinello-Rahman pressure scaling
97 enum class ParrinelloRahmanVelocityScaling
99 No, //!< Do not apply velocity scaling (not a PR-coupling run or step)
100 Diagonal, //!< Apply velocity scaling using a diagonal matrix
101 Full, //!< Apply velocity scaling using a full matrix
102 Count
105 //! Generic callback to the propagator
106 typedef std::function<void(Step)> PropagatorCallback;
107 //! Pointer to generic callback to the propagator
108 typedef std::unique_ptr<PropagatorCallback> PropagatorCallbackPtr;
110 /*! \internal
111 * \brief Propagator element
113 * The propagator element can, through templating, cover the different
114 * propagation types used in NVE MD. The combination of templating, static
115 * functions, and having only the inner-most operations in the static
116 * functions allows to have performance comparable to fused update elements
117 * while keeping easily re-orderable single instructions.
119 * @tparam algorithm The integration types
121 template<IntegrationStep algorithm>
122 class Propagator final : public ISimulatorElement
124 public:
125 //! Constructor
126 Propagator(double timestep,
127 StatePropagatorData* statePropagatorData,
128 const MDAtoms* mdAtoms,
129 gmx_wallcycle* wcycle);
131 /*! \brief Register run function for step / time
133 * @param step The step number
134 * @param time The time
135 * @param registerRunFunction Function allowing to register a run function
137 void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
139 //! No element setup needed
140 void elementSetup() override {}
141 //! No element teardown needed
142 void elementTeardown() override {}
144 //! Set the number of velocity scaling variables
145 void setNumVelocityScalingVariables(int numVelocityScalingVariables);
146 //! Get view on the velocity scaling vector
147 ArrayRef<real> viewOnVelocityScaling();
148 //! Get velocity scaling callback
149 PropagatorCallbackPtr velocityScalingCallback();
151 //! Get view on the full PR scaling matrix
152 ArrayRef<rvec> viewOnPRScalingMatrix();
153 //! Get PR scaling callback
154 PropagatorCallbackPtr prScalingCallback();
156 private:
157 //! The actual propagation
158 template<NumVelocityScalingValues numVelocityScalingValues, ParrinelloRahmanVelocityScaling parrinelloRahmanVelocityScaling>
159 void run();
161 //! The time step
162 const real timestep_;
164 //! Pointer to the micro state
165 StatePropagatorData* statePropagatorData_;
167 //! Whether we're doing single-value velocity scaling
168 bool doSingleVelocityScaling_;
169 //! Wether we're doing group-wise velocity scaling
170 bool doGroupVelocityScaling_;
171 //! The vector of velocity scaling values
172 std::vector<real> velocityScaling_;
173 //! The next velocity scaling step
174 Step scalingStepVelocity_;
176 //! The diagonal of the PR scaling matrix
177 rvec diagPR_;
178 //! The full PR scaling matrix
179 matrix matrixPR_;
180 //! The next PR scaling step
181 Step scalingStepPR_;
183 // Access to ISimulator data
184 //! Atom parameters for this domain.
185 const MDAtoms* mdAtoms_;
186 //! Manages wall cycle accounting.
187 gmx_wallcycle* wcycle_;
190 //! \}
191 } // namespace gmx
193 #endif // GMX_MODULARSIMULATOR_PROPAGATOR_H