Add InteractionDefinitions
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / constrtestdata.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,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.
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 SHAKE and LINCS tests header.
37  *
38  * Contains description and constructor for the test data accumulating object,
39  * declares CPU- and GPU-based functions used to apply SHAKE or LINCS on the
40  * test data.
41  *
42  * \author Artem Zhmurov <zhmurov@gmail.com>
43  * \ingroup module_mdlib
44  */
45
46 #ifndef GMX_MDLIB_TESTS_CONSTRTESTDATA_H
47 #define GMX_MDLIB_TESTS_CONSTRTESTDATA_H
48
49 #include <memory>
50 #include <vector>
51
52 #include "gromacs/gmxlib/nrnb.h"
53 #include "gromacs/gpu_utils/gpu_testutils.h"
54 #include "gromacs/math/paddedvector.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/mdlib/gmx_omp_nthreads.h"
58 #include "gromacs/mdlib/lincs.h"
59 #include "gromacs/mdlib/shake.h"
60 #include "gromacs/mdrunutility/multisim.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/inputrec.h"
63 #include "gromacs/mdtypes/mdatom.h"
64 #include "gromacs/pbcutil/pbc.h"
65 #include "gromacs/topology/idef.h"
66 #include "gromacs/topology/ifunc.h"
67 #include "gromacs/topology/topology.h"
68
69 namespace gmx
70 {
71 namespace test
72 {
73
74 /* \brief
75  * Constraints test data structure.
76  *
77  * Structure to collect all the necessary data, including system coordinates and topology,
78  * constraints information, etc. The structure can be reset and reused.
79  */
80 class ConstraintsTestData
81 {
82 public:
83     //! Human-friendly name for a system
84     std::string title_;
85     //! Number of atoms
86     int numAtoms_;
87     //! Topology
88     gmx_mtop_t mtop_;
89     //! Masses
90     std::vector<real> masses_;
91     //! Inverse masses
92     std::vector<real> invmass_;
93     //! Communication record
94     t_commrec cr_;
95     //! Input record (info that usually in .mdp file)
96     t_inputrec ir_;
97     //! Local topology
98     std::unique_ptr<InteractionDefinitions> idef_;
99     //! MD atoms
100     t_mdatoms md_;
101     //! Multisim data
102     gmx_multisim_t ms_;
103     //! Computational time array (normally used to benchmark performance)
104     t_nrnb nrnb_;
105
106     //! Inverse timestep
107     real invdt_;
108     //! Number of flexible constraints
109     int nflexcon_ = 0;
110     //! Whether the virial should be computed
111     bool computeVirial_;
112     //! Scaled virial
113     tensor virialScaled_;
114     //! Scaled virial (reference values)
115     tensor virialScaledRef_;
116     //! If the free energy is computed
117     bool compute_dHdLambda_;
118     //! For free energy computation
119     real dHdLambda_;
120     //! For free energy computation (reference value)
121     real dHdLambdaRef_;
122
123     //! Coordinates before the timestep
124     PaddedVector<RVec> x_;
125     //! Coordinates after timestep, output for the constraints
126     PaddedVector<RVec> xPrime_;
127     //! Backup for coordinates (for reset)
128     PaddedVector<RVec> xPrime0_;
129     //! Intermediate set of coordinates (normally used for projection correction)
130     PaddedVector<RVec> xPrime2_;
131     //! Velocities
132     PaddedVector<RVec> v_;
133     //! Backup for velocities (for reset)
134     PaddedVector<RVec> v0_;
135
136     //! Constraints data (type1-i1-j1-type2-i2-j2-...)
137     std::vector<int> constraints_;
138     //! Target lengths for all constraint types
139     std::vector<real> constraintsR0_;
140
141     /*! \brief
142      * Constructor for the object with all parameters and variables needed by constraints algorithms.
143      *
144      * This constructor assembles stubs for all the data structures, required to initialize
145      * and apply LINCS and SHAKE constraints. The coordinates and velocities before constraining
146      * are saved to allow for reset. The constraints data are stored for testing after constraints
147      * were applied.
148      *
149      * \param[in]  title                Human-friendly name of the system.
150      * \param[in]  numAtoms             Number of atoms in the system.
151      * \param[in]  masses               Atom masses. Size of this vector should be equal to numAtoms.
152      * \param[in]  constraints          List of constraints, organized in triples of integers.
153      *                                  First integer is the index of type for a constraint, second
154      *                                  and third are the indices of constrained atoms. The types
155      *                                  of constraints should be sequential but not necessarily
156      *                                  start from zero (which is the way they normally are in
157      *                                  GROMACS).
158      * \param[in]  constraintsR0        Target values for bond lengths for bonds of each type. The
159      *                                  size of this vector should be equal to the total number of
160      *                                  unique types in constraints vector.
161      * \param[in]  computeVirial        Whether the virial should be computed.
162      * \param[in]  virialScaledRef      Reference values for scaled virial tensor.
163      * \param[in]  compute_dHdLambda    Whether free energy should be computed.
164      * \param[in]  dHdLambdaRef         Reference value for dHdLambda.
165      * \param[in]  initialTime          Initial time.
166      * \param[in]  timestep             Timestep.
167      * \param[in]  x                    Coordinates before integration step.
168      * \param[in]  xPrime               Coordinates after integration step, but before constraining.
169      * \param[in]  v                    Velocities before constraining.
170      * \param[in]  shakeTolerance       Target tolerance for SHAKE.
171      * \param[in]  shakeUseSOR          Use successive over-relaxation method for SHAKE iterations.
172      *                                  The general formula is:
173      *                                     x_n+1 = (1-omega)*x_n + omega*f(x_n),
174      *                                  where omega = 1 if SOR is off and may be < 1 if SOR is on.
175      * \param[in]  lincsNumIterations   Number of iterations used to compute the inverse matrix.
176      * \param[in]  lincsExpansionOrder  The order for algorithm that adjusts the direction of the
177      *                                  bond after constraints are applied.
178      * \param[in]  lincsWarnAngle       The threshold value for the change in bond angle. When
179      *                                  exceeded the program will issue a warning.
180      *
181      */
182     ConstraintsTestData(const std::string&       title,
183                         int                      numAtoms,
184                         std::vector<real>        masses,
185                         std::vector<int>         constraints,
186                         std::vector<real>        constraintsR0,
187                         bool                     computeVirial,
188                         tensor                   virialScaledRef,
189                         bool                     compute_dHdLambda,
190                         float                    dHdLambdaRef,
191                         real                     initialTime,
192                         real                     timestep,
193                         const std::vector<RVec>& x,
194                         const std::vector<RVec>& xPrime,
195                         const std::vector<RVec>& v,
196                         real                     shakeTolerance,
197                         gmx_bool                 shakeUseSOR,
198                         int                      lincsNumIterations,
199                         int                      lincsExpansionOrder,
200                         real                     lincsWarnAngle);
201
202     /*! \brief
203      * Reset the data structure so it can be reused.
204      *
205      * Set the coordinates and velocities back to their values before
206      * constraining. The scaled virial tensor and dHdLambda are zeroed.
207      *
208      */
209     void reset();
210 };
211
212 } // namespace test
213 } // namespace gmx
214
215 #endif // GMX_MDLIB_TESTS_CONSTRTESTDATA_H