631a98f133782fface5d4f07495e118585b1a6b1
[alexxy/gromacs.git] / src / gromacs / modularsimulator / propagator.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief Declares the propagator element for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  *
41  * This header is only used within the modular simulator module
42  */
43
44 #ifndef GMX_MODULARSIMULATOR_PROPAGATOR_H
45 #define GMX_MODULARSIMULATOR_PROPAGATOR_H
46
47 #include <vector>
48
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/real.h"
52
53 #include "modularsimulatorinterfaces.h"
54
55 struct gmx_wallcycle;
56
57 namespace gmx
58 {
59 class EnergyData;
60 class FreeEnergyPerturbationData;
61 class GlobalCommunicationHelper;
62 class LegacySimulatorData;
63 class MDAtoms;
64 class ModularSimulatorAlgorithmBuilderHelper;
65 class StatePropagatorData;
66
67 //! \addtogroup module_modularsimulator
68 //! \{
69
70 //! Which velocities the thermostat scales
71 enum class ScaleVelocities
72 {
73     PreStepOnly,
74     PreStepAndPostStep
75 };
76
77 //! The different integration types we know about
78 enum class IntegrationStage
79 {
80     PositionsOnly,                        //!< Moves the position vector by the given time step
81     VelocitiesOnly,                       //!< Moves the velocity vector by the given time step
82     LeapFrog,                             //!< Manual fusion of the previous two propagators
83     VelocityVerletPositionsAndVelocities, //!< Manual position (full dt) and velocity (half dt) fusion
84     ScaleVelocities,                      //!< Only scale velocities, don't propagate
85     ScalePositions,                       //!< Only scale positions, don't propagate
86     Count                                 //!< The number of enum entries
87 };
88
89 //! Sets the number of different position scaling values
90 enum class NumPositionScalingValues
91 {
92     None,     //!< No position scaling (either this step or ever)
93     Single,   //!< Single scaling value (either one group or all values =1)
94     Multiple, //!< Multiple scaling values, need to use T-group indices
95     Count     //!< The number of enum entries
96 };
97
98 //! Sets the number of different velocity scaling values
99 enum class NumVelocityScalingValues
100 {
101     None,     //!< No velocity scaling (either this step or ever)
102     Single,   //!< Single T-scaling value (either one group or all values =1)
103     Multiple, //!< Multiple T-scaling values, need to use T-group indices
104     Count
105 };
106
107 //! Sets the type of Parrinello-Rahman pressure scaling
108 enum class ParrinelloRahmanVelocityScaling
109 {
110     No,       //!< Do not apply velocity scaling (not a PR-coupling run or step)
111     Diagonal, //!< Apply velocity scaling using a diagonal matrix
112     Full,     //!< Apply velocity scaling using a full matrix
113     Count
114 };
115
116 /*! \internal
117  * \brief Propagator element
118  *
119  * The propagator element can, through templating, cover the different
120  * propagation types used in NVE MD. The combination of templating, static
121  * functions, and having only the inner-most operations in the static
122  * functions allows to have performance comparable to fused update elements
123  * while keeping easily re-orderable single instructions.
124  *
125  * \tparam integrationStep  The integration types
126  */
127 template<IntegrationStage integrationStage>
128 class Propagator final : public ISimulatorElement
129 {
130 public:
131     //! Constructor
132     Propagator(double               timestep,
133                StatePropagatorData* statePropagatorData,
134                const MDAtoms*       mdAtoms,
135                gmx_wallcycle*       wcycle);
136
137     /*! \brief Register run function for step / time
138      *
139      * This function will determine the required flavor of the run function to be registered
140      * for the current step. In case of the pure scaling integrator stage, it might also skip
141      * the function registration if no scaling is needed.
142      *
143      * \param step                 The step number
144      * \param time                 The time
145      * \param registerRunFunction  Function allowing to register a run function
146      */
147     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
148
149     //! No element setup needed
150     void elementSetup() override {}
151     //! No element teardown needed
152     void elementTeardown() override {}
153
154     //! Set the number of velocity scaling variables
155     void setNumVelocityScalingVariables(int numVelocityScalingVariables, ScaleVelocities scaleVelocities);
156     //! Set the number of position scaling variables
157     void setNumPositionScalingVariables(int numPositionScalingVariables);
158     //! Get view on the scaling vector applied to start of step velocities
159     ArrayRef<real> viewOnStartVelocityScaling();
160     //! Get view on the scaling vector applied to end of step velocities
161     ArrayRef<real> viewOnEndVelocityScaling();
162     //! Get view on the scaling vector applied to the positions
163     ArrayRef<real> viewOnPositionScaling();
164     //! Get velocity scaling callback
165     PropagatorCallback velocityScalingCallback();
166     //! Get position scaling callback
167     PropagatorCallback positionScalingCallback();
168
169     //! Get view on the full PR scaling matrix
170     ArrayRef<rvec> viewOnPRScalingMatrix();
171     //! Get PR scaling callback
172     PropagatorCallback prScalingCallback();
173
174     /*! \brief Factory method implementation
175      *
176      * \param legacySimulatorData  Pointer allowing access to simulator level data
177      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
178      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
179      * \param energyData  Pointer to the \c EnergyData object
180      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
181      * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
182      * \param propagatorTag  The name of the propagator to simplify connection
183      * \param timestep  The time step the propagator uses
184      *
185      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
186      */
187     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
188                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
189                                                     StatePropagatorData*        statePropagatorData,
190                                                     EnergyData*                 energyData,
191                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
192                                                     GlobalCommunicationHelper* globalCommunicationHelper,
193                                                     const PropagatorTag&       propagatorTag,
194                                                     TimeStep                   timestep);
195
196     /*! \brief Factory method implementation
197      *
198      * Version without time step for pure scaling elements
199      *
200      * \param legacySimulatorData  Pointer allowing access to simulator level data
201      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
202      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
203      * \param energyData  Pointer to the \c EnergyData object
204      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
205      * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
206      * \param propagatorTag  The name of the propagator to simplify connection
207      *
208      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
209      */
210     static ISimulatorElement* getElementPointerImpl(LegacySimulatorData* legacySimulatorData,
211                                                     ModularSimulatorAlgorithmBuilderHelper* builderHelper,
212                                                     StatePropagatorData*        statePropagatorData,
213                                                     EnergyData*                 energyData,
214                                                     FreeEnergyPerturbationData* freeEnergyPerturbationData,
215                                                     GlobalCommunicationHelper* globalCommunicationHelper,
216                                                     const PropagatorTag&       propagatorTag);
217
218 private:
219     //! The actual propagation
220     template<NumVelocityScalingValues        numStartVelocityScalingValues,
221              ParrinelloRahmanVelocityScaling parrinelloRahmanVelocityScaling,
222              NumVelocityScalingValues        numEndVelocityScalingValues,
223              NumPositionScalingValues        numPositionScalingValues>
224     void run();
225
226     //! The time step
227     const real timestep_;
228
229     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
230     //! Pointer to the micro state
231     StatePropagatorData* statePropagatorData_;
232
233     //! Whether we're doing single-value velocity scaling (velocities at start of step)
234     bool doSingleStartVelocityScaling_;
235     //! Whether we're doing group-wise velocity scaling (velocities at start of step)
236     bool doGroupStartVelocityScaling_;
237     //! Whether we're doing single-value velocity scaling (velocities at end of step)
238     bool doSingleEndVelocityScaling_;
239     //! Wether we're doing group-wise velocity scaling (velocities at end of step)
240     bool doGroupEndVelocityScaling_;
241     //! Whether we're doing single-value position scaling
242     bool doSinglePositionScaling_;
243     //! Whether we're doing group-wise position scaling
244     bool doGroupPositionScaling_;
245     //! The vector of velocity scaling values
246     std::vector<real> startVelocityScaling_;
247     //! The vector of velocity scaling values
248     std::vector<real> endVelocityScaling_;
249     //! The vector of position scaling values
250     std::vector<real> positionScaling_;
251     //! The next velocity scaling step
252     Step scalingStepVelocity_;
253     //! The next position scaling step
254     Step scalingStepPosition_;
255
256     //! The diagonal of the PR scaling matrix
257     rvec diagPR_;
258     //! The full PR scaling matrix
259     matrix matrixPR_;
260     //! The next PR scaling step
261     Step scalingStepPR_;
262
263     // Access to ISimulator data
264     //! Atom parameters for this domain.
265     const MDAtoms* mdAtoms_;
266     //! Manages wall cycle accounting.
267     gmx_wallcycle* wcycle_;
268 };
269
270 //! \}
271 } // namespace gmx
272
273 #endif // GMX_MODULARSIMULATOR_PROPAGATOR_H