Add InteractionDefinitions
[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,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.
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     int maxType = 0;
149     for (index i = 0; i < ssize(constraints); i += 3)
150     {
151         if (maxType < constraints.at(i))
152         {
153             maxType = constraints.at(i);
154         }
155     }
156     auto& iparams = mtop_.ffparams.iparams;
157     iparams.resize(maxType + 1);
158     for (index i = 0; i < ssize(constraints) / 3; i++)
159     {
160         iparams[constraints.at(3 * i)].constr.dA = constraintsR0.at(constraints.at(3 * i));
161         iparams[constraints.at(3 * i)].constr.dB = constraintsR0.at(constraints.at(3 * i));
162     }
163     idef_ = std::make_unique<InteractionDefinitions>(mtop_.ffparams);
164     for (index i = 0; i < ssize(constraints); i++)
165     {
166         idef_->il[F_CONSTR].iatoms.push_back(constraints.at(i));
167     }
168
169     // Constraints and their parameters (global topology)
170     InteractionList interactionList;
171     interactionList.iatoms.resize(constraints.size());
172     std::copy(constraints.begin(), constraints.end(), interactionList.iatoms.begin());
173     InteractionList interactionListEmpty;
174     interactionListEmpty.iatoms.resize(0);
175
176     gmx_moltype_t molType;
177     molType.atoms.nr             = numAtoms;
178     molType.ilist.at(F_CONSTR)   = interactionList;
179     molType.ilist.at(F_CONSTRNC) = interactionListEmpty;
180     mtop_.moltype.push_back(molType);
181
182     gmx_molblock_t molBlock;
183     molBlock.type = 0;
184     molBlock.nmol = 1;
185     mtop_.molblock.push_back(molBlock);
186
187     mtop_.natoms                      = numAtoms;
188     mtop_.bIntermolecularInteractions = false;
189
190     // Coordinates and velocities
191     x_.resizeWithPadding(numAtoms);
192     xPrime_.resizeWithPadding(numAtoms);
193     xPrime0_.resizeWithPadding(numAtoms);
194     xPrime2_.resizeWithPadding(numAtoms);
195
196     v_.resizeWithPadding(numAtoms);
197     v0_.resizeWithPadding(numAtoms);
198
199     std::copy(x.begin(), x.end(), x_.begin());
200     std::copy(xPrime.begin(), xPrime.end(), xPrime_.begin());
201     std::copy(xPrime.begin(), xPrime.end(), xPrime0_.begin());
202     std::copy(xPrime.begin(), xPrime.end(), xPrime2_.begin());
203
204     std::copy(v.begin(), v.end(), v_.begin());
205     std::copy(v.begin(), v.end(), v0_.begin());
206
207     // SHAKE-specific parameters
208     ir_.shake_tol = shakeTolerance;
209     ir_.bShakeSOR = shakeUseSOR;
210
211     // LINCS-specific parameters
212     ir_.nLincsIter     = lincsNumIterations;
213     ir_.nProjOrder     = lincsExpansionOrder;
214     ir_.LincsWarnAngle = lincsWarnAngle;
215 }
216
217 /*! \brief
218  * Reset the data structure so it can be reused.
219  *
220  * Set the coordinates and velocities back to their values before
221  * constraining. The scaled virial tensor and dHdLambda are zeroed.
222  *
223  */
224 void ConstraintsTestData::reset()
225 {
226     xPrime_  = xPrime0_;
227     xPrime2_ = xPrime0_;
228     v_       = v0_;
229
230     if (computeVirial_)
231     {
232         for (int i = 0; i < DIM; i++)
233         {
234             for (int j = 0; j < DIM; j++)
235             {
236                 virialScaled_[i][j] = 0;
237             }
238         }
239     }
240     dHdLambda_ = 0;
241 }
242
243 } // namespace test
244 } // namespace gmx