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