Use new GPU infrastructure in MDLib tests
[gromacs.git] / src / gromacs / mdlib / tests / settletestdata.cpp
blob158db19fe7442c0323d86e421fc0767405d20f03
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief Defines the class that accumulates SETTLE test data.
38 * \author Mark Abraham <mark.j.abraham@gmail.com>
39 * \author Artem Zhmurov <zhmurov@gmail.com>
41 * \ingroup module_mdlib
43 #include "gmxpre.h"
45 #include "settletestdata.h"
47 #include <algorithm>
48 #include <unordered_map>
49 #include <vector>
51 #include <gtest/gtest.h>
53 #include "gromacs/gpu_utils/gpu_utils.h"
54 #include "gromacs/math/paddedvector.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/mdlib/settle.h"
58 #include "gromacs/mdtypes/mdatom.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/topology/idef.h"
61 #include "gromacs/topology/ifunc.h"
62 #include "gromacs/topology/topology.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/stringutil.h"
65 #include "gromacs/utility/unique_cptr.h"
67 #include "gromacs/mdlib/tests/watersystem.h"
68 #include "testutils/testasserts.h"
70 namespace gmx
72 namespace test
75 SettleTestData::SettleTestData(int numSettles) :
76 numSettles_(numSettles),
77 x_(c_waterPositions.size()),
78 xPrime_(c_waterPositions.size()),
79 v_(c_waterPositions.size())
81 // Initialize coordinates and velocities from the constant set of coordinates
82 std::copy(c_waterPositions.begin(), c_waterPositions.end(), x_.begin());
83 std::copy(c_waterPositions.begin(), c_waterPositions.end(), xPrime_.begin());
85 // Perturb the atom positions, to appear like an
86 // "update," and where there is definitely constraining
87 // work to do.
88 const real deltas[] = { 0.01, -0.01, +0.02, -0.02 };
89 int i = 0;
90 for (auto& xPrime : xPrime_)
92 xPrime[XX] += deltas[i % 4];
93 ++i;
94 xPrime[YY] += deltas[i % 4];
95 ++i;
96 xPrime[ZZ] += deltas[i % 4];
97 ++i;
99 std::fill(v_.begin(), v_.end(), RVec{ 0.0, 0.0, 0.0 });
101 // Set up the topology.
102 const int settleType = 0;
103 mtop_.moltype.resize(1);
104 mtop_.molblock.resize(1);
105 mtop_.molblock[0].type = 0;
106 std::vector<int>& iatoms = mtop_.moltype[0].ilist[F_SETTLE].iatoms;
107 for (int i = 0; i < numSettles; ++i)
109 iatoms.push_back(settleType);
110 iatoms.push_back(i * atomsPerSettle_ + 0);
111 iatoms.push_back(i * atomsPerSettle_ + 1);
112 iatoms.push_back(i * atomsPerSettle_ + 2);
115 // Set up the SETTLE parameters.
116 t_iparams iparams;
117 iparams.settle.doh = dOH_;
118 iparams.settle.dhh = dHH_;
119 mtop_.ffparams.iparams.push_back(iparams);
121 // Set up the masses.
122 mtop_.moltype[0].atoms.atom =
123 static_cast<t_atom*>(calloc(numSettles * atomsPerSettle_, sizeof(t_atom)));
124 numAtoms_ = numSettles * atomsPerSettle_;
125 masses_.resize(numAtoms_);
126 inverseMasses_.resize(numAtoms_);
127 for (int i = 0; i < numSettles; ++i)
129 masses_[i * atomsPerSettle_ + 0] = oxygenMass_;
130 masses_[i * atomsPerSettle_ + 1] = hydrogenMass_;
131 masses_[i * atomsPerSettle_ + 2] = hydrogenMass_;
133 inverseMasses_[i * atomsPerSettle_ + 0] = 1.0 / oxygenMass_;
134 inverseMasses_[i * atomsPerSettle_ + 1] = 1.0 / hydrogenMass_;
135 inverseMasses_[i * atomsPerSettle_ + 2] = 1.0 / hydrogenMass_;
137 mtop_.moltype[0].atoms.atom[i * atomsPerSettle_ + 0].m = oxygenMass_;
138 mtop_.moltype[0].atoms.atom[i * atomsPerSettle_ + 1].m = hydrogenMass_;
139 mtop_.moltype[0].atoms.atom[i * atomsPerSettle_ + 2].m = hydrogenMass_;
142 idef_ = std::make_unique<InteractionDefinitions>(mtop_.ffparams);
143 idef_->il[F_SETTLE] = mtop_.moltype[0].ilist[F_SETTLE];
146 SettleTestData::~SettleTestData() {}
148 } // namespace test
149 } // namespace gmx