Merge branch 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) :
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 {
82     mdAtoms_.nr = numAtoms_;
83
84     for (int i = 0; i < numAtoms_; i++)
85     {
86         // Typical PBC box size is tens of nanometers
87         x_[i][XX] = (i%21)*1.0;
88         x_[i][YY] = 6.5 + (i%13)*(-1.0);
89         x_[i][ZZ] = (i%32)*(0.0);
90
91         for (int d = 0; d < DIM; d++)
92         {
93             xPrime_[i][d] = 0.0;
94             // Thermal velocity is ~1 nm/ps (|v0| = 1-2 nm/ps)
95             v_[i][d] = v0[d];
96             // TODO Check what value typical MD forces have (now ~ 1 kJ/mol/nm)
97             f_[i][d] = f0[d];
98
99             x0_[i][d] = x_[i][d];
100             v0_[i][d] = v_[i][d];
101         }
102         // Atom masses are ~1-100 g/mol
103         inverseMasses_[i] = 1.0/(1.0 + i%100);
104         for (int d = 0; d < DIM; d++)
105         {
106             inverseMassesPerDim_[i][d] = inverseMasses_[i];
107         }
108     }
109     mdAtoms_.invmass       = inverseMasses_.data();
110     mdAtoms_.invMassPerDim = as_rvec_array(inverseMassesPerDim_.data());
111
112     // Data needed for current CPU-based implementation
113     inputRecord_.eI      = eiMD;
114     inputRecord_.delta_t = timestep_;
115     inputRecord_.etc     = etcNO;
116     inputRecord_.epc     = epcNO;
117
118     state_.flags = 0;
119
120     state_.box[XX][XX] = 10.0;
121     state_.box[XX][YY] = 0.0;
122     state_.box[XX][ZZ] = 0.0;
123
124     state_.box[YY][XX] = 0.0;
125     state_.box[YY][YY] = 10.0;
126     state_.box[YY][ZZ] = 0.0;
127
128     state_.box[ZZ][XX] = 0.0;
129     state_.box[ZZ][YY] = 0.0;
130     state_.box[ZZ][ZZ] = 10.0;
131
132     kineticEnergyData_.bNEMD            = false;
133     kineticEnergyData_.cosacc.cos_accel = 0.0;
134     t_grp_tcstat temperatureCouplingGroupData;
135     temperatureCouplingGroupData.lambda = 1;
136     kineticEnergyData_.tcstat.emplace_back(temperatureCouplingGroupData);
137
138     kineticEnergyData_.nthreads = 1;
139     snew(kineticEnergyData_.ekin_work_alloc, kineticEnergyData_.nthreads);
140     snew(kineticEnergyData_.ekin_work, kineticEnergyData_.nthreads);
141     snew(kineticEnergyData_.dekindl_work, kineticEnergyData_.nthreads);
142
143     mdAtoms_.homenr                   = numAtoms_;
144     mdAtoms_.haveVsites               = false;
145     mdAtoms_.havePartiallyFrozenAtoms = false;
146     snew(mdAtoms_.cTC, numAtoms_);
147     for (int i = 0; i < numAtoms_; i++)
148     {
149         mdAtoms_.cTC[i] = 0;
150     }
151
152     prVScalingMatrix_[XX][XX] = 1.0;
153     prVScalingMatrix_[XX][YY] = 0.0;
154     prVScalingMatrix_[XX][ZZ] = 0.0;
155
156     prVScalingMatrix_[YY][XX] = 0.0;
157     prVScalingMatrix_[YY][YY] = 1.0;
158     prVScalingMatrix_[YY][ZZ] = 0.0;
159
160     prVScalingMatrix_[ZZ][XX] = 0.0;
161     prVScalingMatrix_[ZZ][YY] = 0.0;
162     prVScalingMatrix_[ZZ][ZZ] = 1.0;
163
164     update_ = std::make_unique<Update>(&inputRecord_, nullptr);
165     update_->setNumAtoms(numAtoms);
166 }
167
168 LeapFrogTestData::~LeapFrogTestData()
169 {
170     sfree(mdAtoms_.cTC);
171 }
172
173 }  // namespace test
174 }  // namespace gmx