Apply clang-format to source tree
[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, 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 /*! \libinternal \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
42 #ifndef GMX_MODULARSIMULATOR_PROPAGATOR_H
43 #define GMX_MODULARSIMULATOR_PROPAGATOR_H
44
45 #include <vector>
46
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/real.h"
50
51 #include "modularsimulatorinterfaces.h"
52
53 struct gmx_wallcycle;
54
55 namespace gmx
56 {
57 class MDAtoms;
58 class StatePropagatorData;
59
60 //! \addtogroup module_modularsimulator
61 //! \{
62
63 /*! \brief The different integration types we know about
64  *
65  * PositionsOnly:
66  *   Moves the position vector by the given time step
67  * VelocitiesOnly:
68  *   Moves the velocity vector by the given time step
69  * LeapFrog:
70  *   Is a manual fusion of the previous two propagators
71  * VelocityVerletPositionsAndVelocities:
72  *   Is a manual fusion of VelocitiesOnly and PositionsOnly,
73  *   where VelocitiesOnly is only propagated by half the
74  *   time step of the positions.
75  */
76 enum class IntegrationStep
77 {
78     PositionsOnly,
79     VelocitiesOnly,
80     LeapFrog,
81     VelocityVerletPositionsAndVelocities,
82     Count
83 };
84
85 //! Sets the number of different velocity scaling values
86 enum class NumVelocityScalingValues
87 {
88     None,     //!< No velocity scaling (either this step or ever)
89     Single,   //!< Single T-scaling value (either one group or all values =1)
90     Multiple, //!< Multiple T-scaling values, need to use T-group indices
91     Count
92 };
93
94 //! Sets the type of Parrinello-Rahman pressure scaling
95 enum class ParrinelloRahmanVelocityScaling
96 {
97     No,       //!< Do not apply velocity scaling (not a PR-coupling run or step)
98     Diagonal, //!< Apply velocity scaling using a diagonal matrix
99     Full,     //!< Apply velocity scaling using a full matrix
100     Count
101 };
102
103 //! Generic callback to the propagator
104 typedef std::function<void(Step)> PropagatorCallback;
105 //! Pointer to generic callback to the propagator
106 typedef std::unique_ptr<PropagatorCallback> PropagatorCallbackPtr;
107
108 /*! \libinternal
109  * \brief Propagator element
110  *
111  * The propagator element can, through templating, cover the different
112  * propagation types used in NVE MD. The combination of templating, static
113  * functions, and having only the inner-most operations in the static
114  * functions allows to have performance comparable to fused update elements
115  * while keeping easily re-orderable single instructions.
116  *
117  * \todo: Get rid of updateVelocities2() once we don't require identical
118  *        reproduction of do_md() results.
119  *
120  * @tparam algorithm  The integration types
121  */
122 template<IntegrationStep algorithm>
123 class Propagator final : public ISimulatorElement
124 {
125 public:
126     //! Constructor
127     Propagator(double               timestep,
128                StatePropagatorData* statePropagatorData,
129                const MDAtoms*       mdAtoms,
130                gmx_wallcycle*       wcycle);
131
132     /*! \brief Register run function for step / time
133      *
134      * @param step                 The step number
135      * @param time                 The time
136      * @param registerRunFunction  Function allowing to register a run function
137      */
138     void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
139
140     //! No element setup needed
141     void elementSetup() override {}
142     //! No element teardown needed
143     void elementTeardown() override {}
144
145     //! Set the number of velocity scaling variables
146     void setNumVelocityScalingVariables(int numVelocityScalingVariables);
147     //! Get view on the velocity scaling vector
148     ArrayRef<real> viewOnVelocityScaling();
149     //! Get velocity scaling callback
150     PropagatorCallbackPtr velocityScalingCallback();
151
152     //! Get view on the full PR scaling matrix
153     ArrayRef<rvec> viewOnPRScalingMatrix();
154     //! Get PR scaling callback
155     PropagatorCallbackPtr prScalingCallback();
156
157 private:
158     //! The actual propagation
159     template<NumVelocityScalingValues numVelocityScalingValues, ParrinelloRahmanVelocityScaling parrinelloRahmanVelocityScaling>
160     void run();
161
162     //! The time step
163     const real timestep_;
164
165     //! Pointer to the micro state
166     StatePropagatorData* statePropagatorData_;
167
168     //! Whether we're doing single-value velocity scaling
169     bool doSingleVelocityScaling;
170     //! Wether we're doing group-wise velocity scaling
171     bool doGroupVelocityScaling;
172     //! The vector of velocity scaling values
173     std::vector<real> velocityScaling_;
174     //! The next velocity scaling step
175     Step scalingStepVelocity_;
176
177     //! The diagonal of the PR scaling matrix
178     rvec diagPR;
179     //! The full PR scaling matrix
180     matrix matrixPR;
181     //! The next PR scaling step
182     Step scalingStepPR_;
183
184     // Access to ISimulator data
185     //! Atom parameters for this domain.
186     const MDAtoms* mdAtoms_;
187     //! Manages wall cycle accounting.
188     gmx_wallcycle* wcycle_;
189 };
190
191 //! \}
192 } // namespace gmx
193
194 #endif // GMX_MODULARSIMULATOR_PROPAGATOR_H