Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / constrtestdata.cpp
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.
37  *
38  * \todo Better tests for virial are needed.
39  * \todo Tests for bigger systems to test threads synchronization,
40  *       reduction, etc. on the GPU.
41  * \todo Tests for algorithms for derivatives.
42  * \todo Free-energy perturbation tests
43  *
44  * \author Artem Zhmurov <zhmurov@gmail.com>
45  * \ingroup module_mdlib
46  */
47
48 #include "gmxpre.h"
49
50 #include "constrtestdata.h"
51
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/utility/stringutil.h"
54
55 namespace gmx
56 {
57 namespace test
58 {
59
60 ConstraintsTestData::ConstraintsTestData(const std::string&       title,
61                                          int                      numAtoms,
62                                          std::vector<real>        masses,
63                                          std::vector<int>         constraints,
64                                          std::vector<real>        constraintsR0,
65                                          bool                     computeVirial,
66                                          tensor                   virialScaledRef,
67                                          bool                     compute_dHdLambda,
68                                          float                    dHdLambdaRef,
69                                          real                     initialTime,
70                                          real                     timestep,
71                                          const std::vector<RVec>& x,
72                                          const std::vector<RVec>& xPrime,
73                                          const std::vector<RVec>& v,
74                                          real                     shakeTolerance,
75                                          gmx_bool                 shakeUseSOR,
76                                          int                      lincsNumIterations,
77                                          int                      lincsExpansionOrder,
78                                          real                     lincsWarnAngle)
79 {
80     title_    = title;    // Human-friendly name of the system
81     numAtoms_ = numAtoms; // Number of atoms
82
83     // Masses of atoms
84     masses_ = masses;
85     invmass_.resize(numAtoms); // Vector of inverse masses
86
87     for (int i = 0; i < numAtoms; i++)
88     {
89         invmass_[i] = 1.0 / masses.at(i);
90     }
91
92     // Saving constraints to check if they are satisfied after algorithm was applied
93     constraints_   = constraints;   // Constraints indices (in type-i-j format)
94     constraintsR0_ = constraintsR0; // Equilibrium distances for each type of constraint
95
96     invdt_ = 1.0 / timestep; // Inverse timestep
97
98     // Communication record
99     cr_.nnodes = 1;
100     cr_.dd     = nullptr;
101
102     // Multisim data
103     ms_.sim  = 0;
104     ms_.nsim = 1;
105
106     // Input record - data that usually comes from configuration file (.mdp)
107     ir_.efep    = 0;
108     ir_.init_t  = initialTime;
109     ir_.delta_t = timestep;
110     ir_.eI      = 0;
111
112     // MD atoms data
113     md_.nMassPerturbed = 0;
114     md_.lambda         = 0.0;
115     md_.invmass        = invmass_.data();
116     md_.nr             = numAtoms;
117     md_.homenr         = numAtoms;
118
119     // Virial evaluation
120     computeVirial_ = computeVirial;
121     if (computeVirial)
122     {
123         for (int i = 0; i < DIM; i++)
124         {
125             for (int j = 0; j < DIM; j++)
126             {
127                 virialScaled_[i][j]    = 0;
128                 virialScaledRef_[i][j] = virialScaledRef[i][j];
129             }
130         }
131     }
132
133
134     // Free energy evaluation
135     compute_dHdLambda_ = compute_dHdLambda;
136     dHdLambda_         = 0;
137     if (compute_dHdLambda_)
138     {
139         ir_.efep      = efepYES;
140         dHdLambdaRef_ = dHdLambdaRef;
141     }
142     else
143     {
144         ir_.efep      = efepNO;
145         dHdLambdaRef_ = 0;
146     }
147
148     // Constraints and their parameters (local topology)
149     for (int i = 0; i < F_NRE; i++)
150     {
151         idef_.il[i].nr = 0;
152     }
153     idef_.il[F_CONSTR].nr = constraints.size();
154
155     snew(idef_.il[F_CONSTR].iatoms, constraints.size());
156     int maxType = 0;
157     for (index i = 0; i < ssize(constraints); i++)
158     {
159         if (i % 3 == 0)
160         {
161             if (maxType < constraints.at(i))
162             {
163                 maxType = constraints.at(i);
164             }
165         }
166         idef_.il[F_CONSTR].iatoms[i] = constraints.at(i);
167     }
168     snew(idef_.iparams, maxType + 1);
169     for (index i = 0; i < ssize(constraints) / 3; i++)
170     {
171         idef_.iparams[constraints.at(3 * i)].constr.dA = constraintsR0.at(constraints.at(3 * i));
172         idef_.iparams[constraints.at(3 * i)].constr.dB = constraintsR0.at(constraints.at(3 * i));
173     }
174
175     // Constraints and their parameters (global topology)
176     InteractionList interactionList;
177     interactionList.iatoms.resize(constraints.size());
178     std::copy(constraints.begin(), constraints.end(), interactionList.iatoms.begin());
179     InteractionList interactionListEmpty;
180     interactionListEmpty.iatoms.resize(0);
181
182     gmx_moltype_t molType;
183     molType.atoms.nr             = numAtoms;
184     molType.ilist.at(F_CONSTR)   = interactionList;
185     molType.ilist.at(F_CONSTRNC) = interactionListEmpty;
186     mtop_.moltype.push_back(molType);
187
188     gmx_molblock_t molBlock;
189     molBlock.type = 0;
190     molBlock.nmol = 1;
191     mtop_.molblock.push_back(molBlock);
192
193     mtop_.natoms = numAtoms;
194     mtop_.ffparams.iparams.resize(maxType + 1);
195     for (int i = 0; i <= maxType; i++)
196     {
197         mtop_.ffparams.iparams.at(i) = idef_.iparams[i];
198     }
199     mtop_.bIntermolecularInteractions = false;
200
201     // Coordinates and velocities
202     x_.resizeWithPadding(numAtoms);
203     xPrime_.resizeWithPadding(numAtoms);
204     xPrime0_.resizeWithPadding(numAtoms);
205     xPrime2_.resizeWithPadding(numAtoms);
206
207     v_.resizeWithPadding(numAtoms);
208     v0_.resizeWithPadding(numAtoms);
209
210     std::copy(x.begin(), x.end(), x_.begin());
211     std::copy(xPrime.begin(), xPrime.end(), xPrime_.begin());
212     std::copy(xPrime.begin(), xPrime.end(), xPrime0_.begin());
213     std::copy(xPrime.begin(), xPrime.end(), xPrime2_.begin());
214
215     std::copy(v.begin(), v.end(), v_.begin());
216     std::copy(v.begin(), v.end(), v0_.begin());
217
218     // SHAKE-specific parameters
219     ir_.shake_tol = shakeTolerance;
220     ir_.bShakeSOR = shakeUseSOR;
221
222     // LINCS-specific parameters
223     ir_.nLincsIter     = lincsNumIterations;
224     ir_.nProjOrder     = lincsExpansionOrder;
225     ir_.LincsWarnAngle = lincsWarnAngle;
226 }
227
228 /*! \brief
229  * Reset the data structure so it can be reused.
230  *
231  * Set the coordinates and velocities back to their values before
232  * constraining. The scaled virial tensor and dHdLambda are zeroed.
233  *
234  */
235 void ConstraintsTestData::reset()
236 {
237     xPrime_  = xPrime0_;
238     xPrime2_ = xPrime0_;
239     v_       = v0_;
240
241     if (computeVirial_)
242     {
243         for (int i = 0; i < DIM; i++)
244         {
245             for (int j = 0; j < DIM; j++)
246             {
247                 virialScaled_[i][j] = 0;
248             }
249         }
250     }
251     dHdLambda_ = 0;
252 }
253
254 /*! \brief
255  * Cleaning up the memory.
256  */
257 ConstraintsTestData::~ConstraintsTestData()
258 {
259     sfree(idef_.il[F_CONSTR].iatoms);
260     sfree(idef_.iparams);
261 }
262
263 } // namespace test
264 } // namespace gmx