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