3f66855450e8cdec66b8e57c54b1db109c8429ac
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / leapfrogtestdata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020,2021, 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 Defines the class to accumulate the data needed for the Leap-Frog integrator tests
37  *
38  * \author Artem Zhmurov <zhmurov@gmail.com>
39  * \ingroup module_mdlib
40  */
41 #include "gmxpre.h"
42
43 #include "leapfrogtestdata.h"
44
45 #include <assert.h>
46
47 #include <cmath>
48
49 #include <algorithm>
50 #include <unordered_map>
51 #include <vector>
52
53 #include <gtest/gtest.h>
54
55 #include "gromacs/gpu_utils/gpu_utils.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/utility/smalloc.h"
60 #include "gromacs/utility/stringutil.h"
61
62 #include "testutils/refdata.h"
63 #include "testutils/testasserts.h"
64
65 namespace gmx
66 {
67 namespace test
68 {
69
70 LeapFrogTestData::LeapFrogTestData(int        numAtoms,
71                                    real       timestep,
72                                    const rvec v0,
73                                    const rvec f0,
74                                    int        numTCoupleGroups,
75                                    int        nstpcouple) :
76     numAtoms_(numAtoms),
77     timestep_(timestep),
78     x0_(numAtoms),
79     x_(numAtoms),
80     xPrime_(numAtoms),
81     v0_(numAtoms),
82     v_(numAtoms),
83     f_(numAtoms),
84     inverseMasses_(numAtoms),
85     inverseMassesPerDim_(numAtoms),
86     kineticEnergyData_(numTCoupleGroups == 0 ? 1 : numTCoupleGroups, 0.0, 1),
87     numTCoupleGroups_(numTCoupleGroups)
88 {
89     mdAtoms_.nr = numAtoms_;
90
91     for (int i = 0; i < numAtoms_; i++)
92     {
93         // Typical PBC box size is tens of nanometers
94         x_[i][XX] = (i % 21) * 1.0;
95         x_[i][YY] = 6.5 + (i % 13) * (-1.0);
96         x_[i][ZZ] = (i % 32) * (0.0);
97
98         for (int d = 0; d < DIM; d++)
99         {
100             xPrime_[i][d] = 0.0;
101             // Thermal velocity is ~1 nm/ps (|v0| = 1-2 nm/ps)
102             v_[i][d] = v0[d];
103             // TODO Check what value typical MD forces have (now ~ 1 kJ/mol/nm)
104             f_[i][d] = f0[d];
105
106             x0_[i][d] = x_[i][d];
107             v0_[i][d] = v_[i][d];
108         }
109         // Atom masses are ~1-100 g/mol
110         inverseMasses_[i] = 1.0 / (1.0 + i % 100);
111         for (int d = 0; d < DIM; d++)
112         {
113             inverseMassesPerDim_[i][d] = inverseMasses_[i];
114         }
115     }
116     mdAtoms_.invmass       = inverseMasses_.data();
117     mdAtoms_.invMassPerDim = as_rvec_array(inverseMassesPerDim_.data());
118
119     // Temperature coupling
120     snew(mdAtoms_.cTC, numAtoms_);
121
122     // To do temperature coupling at each step
123     inputRecord_.nsttcouple = 1;
124
125     if (numTCoupleGroups_ == 0)
126     {
127         inputRecord_.etc = TemperatureCoupling::No;
128         for (int i = 0; i < numAtoms_; i++)
129         {
130             mdAtoms_.cTC[i] = 0;
131         }
132         t_grp_tcstat temperatureCouplingGroupData;
133         temperatureCouplingGroupData.lambda = 1.0;
134         kineticEnergyData_.tcstat[0]        = temperatureCouplingGroupData;
135     }
136     else
137     {
138         inputRecord_.etc = TemperatureCoupling::Yes;
139         for (int i = 0; i < numAtoms_; i++)
140         {
141             mdAtoms_.cTC[i] = i % numTCoupleGroups_;
142         }
143         for (int i = 0; i < numTCoupleGroups_; i++)
144         {
145             real         tCoupleLambda = 1.0 - (i + 1.0) / 10.0;
146             t_grp_tcstat temperatureCouplingGroupData;
147             temperatureCouplingGroupData.lambda = tCoupleLambda;
148             kineticEnergyData_.tcstat[i]        = temperatureCouplingGroupData;
149         }
150     }
151
152     inputRecord_.eI      = IntegrationAlgorithm::MD;
153     inputRecord_.delta_t = timestep_;
154
155     state_.flags = 0;
156
157     state_.box[XX][XX] = 10.0;
158     state_.box[XX][YY] = 0.0;
159     state_.box[XX][ZZ] = 0.0;
160
161     state_.box[YY][XX] = 0.0;
162     state_.box[YY][YY] = 10.0;
163     state_.box[YY][ZZ] = 0.0;
164
165     state_.box[ZZ][XX] = 0.0;
166     state_.box[ZZ][YY] = 0.0;
167     state_.box[ZZ][ZZ] = 10.0;
168
169     mdAtoms_.homenr                   = numAtoms_;
170     mdAtoms_.haveVsites               = false;
171     mdAtoms_.havePartiallyFrozenAtoms = false;
172     mdAtoms_.cFREEZE                  = nullptr;
173     mdAtoms_.ptype                    = nullptr;
174
175     update_ = std::make_unique<Update>(inputRecord_, nullptr);
176     update_->updateAfterPartition(numAtoms,
177                                   gmx::ArrayRef<const unsigned short>(),
178                                   gmx::arrayRefFromArray(mdAtoms_.cTC, mdAtoms_.nr));
179
180     doPressureCouple_ = (nstpcouple != 0);
181
182     if (doPressureCouple_)
183     {
184         inputRecord_.epc        = PressureCoupling::ParrinelloRahman;
185         inputRecord_.nstpcouple = nstpcouple;
186         dtPressureCouple_       = inputRecord_.nstpcouple * inputRecord_.delta_t;
187
188         velocityScalingMatrix_[XX][XX] = 1.2;
189         velocityScalingMatrix_[XX][YY] = 0.0;
190         velocityScalingMatrix_[XX][ZZ] = 0.0;
191
192         velocityScalingMatrix_[YY][XX] = 0.0;
193         velocityScalingMatrix_[YY][YY] = 0.8;
194         velocityScalingMatrix_[YY][ZZ] = 0.0;
195
196         velocityScalingMatrix_[ZZ][XX] = 0.0;
197         velocityScalingMatrix_[ZZ][YY] = 0.0;
198         velocityScalingMatrix_[ZZ][ZZ] = 0.9;
199     }
200     else
201     {
202         inputRecord_.epc               = PressureCoupling::No;
203         velocityScalingMatrix_[XX][XX] = 1.0;
204         velocityScalingMatrix_[XX][YY] = 0.0;
205         velocityScalingMatrix_[XX][ZZ] = 0.0;
206
207         velocityScalingMatrix_[YY][XX] = 0.0;
208         velocityScalingMatrix_[YY][YY] = 1.0;
209         velocityScalingMatrix_[YY][ZZ] = 0.0;
210
211         velocityScalingMatrix_[ZZ][XX] = 0.0;
212         velocityScalingMatrix_[ZZ][YY] = 0.0;
213         velocityScalingMatrix_[ZZ][ZZ] = 1.0;
214     }
215 }
216
217 LeapFrogTestData::~LeapFrogTestData()
218 {
219     sfree(mdAtoms_.cTC);
220 }
221
222 } // namespace test
223 } // namespace gmx