Merge remote-tracking branch 'origin/release-2019'
[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, 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, real timestep, const rvec v0, const rvec f0, int numTCoupleGroups, int nstpcouple) :
71     numAtoms_(numAtoms),
72     timestep_(timestep),
73     x0_(numAtoms),
74     x_(numAtoms),
75     xPrime_(numAtoms),
76     v0_(numAtoms),
77     v_(numAtoms),
78     f_(numAtoms),
79     inverseMasses_(numAtoms),
80     inverseMassesPerDim_(numAtoms),
81     numTCoupleGroups_(numTCoupleGroups)
82 {
83     mdAtoms_.nr = numAtoms_;
84
85     for (int i = 0; i < numAtoms_; i++)
86     {
87         // Typical PBC box size is tens of nanometers
88         x_[i][XX] = (i%21)*1.0;
89         x_[i][YY] = 6.5 + (i%13)*(-1.0);
90         x_[i][ZZ] = (i%32)*(0.0);
91
92         for (int d = 0; d < DIM; d++)
93         {
94             xPrime_[i][d] = 0.0;
95             // Thermal velocity is ~1 nm/ps (|v0| = 1-2 nm/ps)
96             v_[i][d] = v0[d];
97             // TODO Check what value typical MD forces have (now ~ 1 kJ/mol/nm)
98             f_[i][d] = f0[d];
99
100             x0_[i][d] = x_[i][d];
101             v0_[i][d] = v_[i][d];
102         }
103         // Atom masses are ~1-100 g/mol
104         inverseMasses_[i] = 1.0/(1.0 + i%100);
105         for (int d = 0; d < DIM; d++)
106         {
107             inverseMassesPerDim_[i][d] = inverseMasses_[i];
108         }
109     }
110     mdAtoms_.invmass       = inverseMasses_.data();
111     mdAtoms_.invMassPerDim = as_rvec_array(inverseMassesPerDim_.data());
112
113     // Temperature coupling
114     snew(mdAtoms_.cTC, numAtoms_);
115
116     // To do temperature coupling at each step
117     inputRecord_.nsttcouple = 1;
118
119     if (numTCoupleGroups_ == 0)
120     {
121         inputRecord_.etc = etcNO;
122         for (int i = 0; i < numAtoms_; i++)
123         {
124             mdAtoms_.cTC[i] = 0;
125         }
126         kineticEnergyData_.ngtc = 1;
127         t_grp_tcstat temperatureCouplingGroupData;
128         temperatureCouplingGroupData.lambda = 1.0;
129         kineticEnergyData_.tcstat.emplace_back(temperatureCouplingGroupData);
130     }
131     else
132     {
133         inputRecord_.etc = etcYES;
134         for (int i = 0; i < numAtoms_; i++)
135         {
136             mdAtoms_.cTC[i] = i % numTCoupleGroups_;
137         }
138         kineticEnergyData_.ngtc = numTCoupleGroups_;
139         for (int i = 0; i < numTCoupleGroups; i++)
140         {
141             real         tCoupleLambda = 1.0 - (i + 1.0)/10.0;
142             t_grp_tcstat temperatureCouplingGroupData;
143             temperatureCouplingGroupData.lambda = tCoupleLambda;
144             kineticEnergyData_.tcstat.emplace_back(temperatureCouplingGroupData);
145         }
146     }
147
148     inputRecord_.eI      = eiMD;
149     inputRecord_.delta_t = timestep_;
150
151     state_.flags = 0;
152
153     state_.box[XX][XX] = 10.0;
154     state_.box[XX][YY] = 0.0;
155     state_.box[XX][ZZ] = 0.0;
156
157     state_.box[YY][XX] = 0.0;
158     state_.box[YY][YY] = 10.0;
159     state_.box[YY][ZZ] = 0.0;
160
161     state_.box[ZZ][XX] = 0.0;
162     state_.box[ZZ][YY] = 0.0;
163     state_.box[ZZ][ZZ] = 10.0;
164
165     kineticEnergyData_.bNEMD            = false;
166     kineticEnergyData_.cosacc.cos_accel = 0.0;
167
168     kineticEnergyData_.nthreads = 1;
169     snew(kineticEnergyData_.ekin_work_alloc, kineticEnergyData_.nthreads);
170     snew(kineticEnergyData_.ekin_work, kineticEnergyData_.nthreads);
171     snew(kineticEnergyData_.dekindl_work, kineticEnergyData_.nthreads);
172
173     mdAtoms_.homenr                   = numAtoms_;
174     mdAtoms_.haveVsites               = false;
175     mdAtoms_.havePartiallyFrozenAtoms = false;
176     mdAtoms_.cFREEZE                  = nullptr;
177
178     update_ = std::make_unique<Update>(&inputRecord_, nullptr);
179     update_->setNumAtoms(numAtoms);
180
181     doPressureCouple_ = (nstpcouple != 0);
182
183     if (doPressureCouple_)
184     {
185         inputRecord_.epc        = epcPARRINELLORAHMAN;
186         inputRecord_.nstpcouple = nstpcouple;
187         dtPressureCouple_       = inputRecord_.nstpcouple*inputRecord_.delta_t;
188
189         velocityScalingMatrix_[XX][XX] = 1.2;
190         velocityScalingMatrix_[XX][YY] = 0.0;
191         velocityScalingMatrix_[XX][ZZ] = 0.0;
192
193         velocityScalingMatrix_[YY][XX] = 0.0;
194         velocityScalingMatrix_[YY][YY] = 0.8;
195         velocityScalingMatrix_[YY][ZZ] = 0.0;
196
197         velocityScalingMatrix_[ZZ][XX] = 0.0;
198         velocityScalingMatrix_[ZZ][YY] = 0.0;
199         velocityScalingMatrix_[ZZ][ZZ] = 0.9;
200     }
201     else
202     {
203         inputRecord_.epc               = epcNO;
204         velocityScalingMatrix_[XX][XX] = 1.0;
205         velocityScalingMatrix_[XX][YY] = 0.0;
206         velocityScalingMatrix_[XX][ZZ] = 0.0;
207
208         velocityScalingMatrix_[YY][XX] = 0.0;
209         velocityScalingMatrix_[YY][YY] = 1.0;
210         velocityScalingMatrix_[YY][ZZ] = 0.0;
211
212         velocityScalingMatrix_[ZZ][XX] = 0.0;
213         velocityScalingMatrix_[ZZ][YY] = 0.0;
214         velocityScalingMatrix_[ZZ][ZZ] = 1.0;
215     }
216
217 }
218
219 LeapFrogTestData::~LeapFrogTestData()
220 {
221     sfree(mdAtoms_.cTC);
222 }
223
224 }  // namespace test
225 }  // namespace gmx