2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2018,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.
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.
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
45 #include "settletestdata.h"
48 #include <unordered_map>
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"
75 SettleTestData::SettleTestData(int numSettles) :
76 x_(c_waterPositions.size()),
77 xPrime_(c_waterPositions.size()),
78 v_(c_waterPositions.size())
80 // Initialize coordinates and velocities from the constant set of coordinates
81 std::copy(c_waterPositions.begin(), c_waterPositions.end(), x_.begin());
82 std::copy(c_waterPositions.begin(), c_waterPositions.end(), xPrime_.begin());
84 // Perturb the atom positions, to appear like an
85 // "update," and where there is definitely constraining
87 const real deltas[] = { 0.01, -0.01, +0.02, -0.02 };
89 for (auto &xPrime : xPrime_)
91 xPrime[XX] += deltas[i % 4];
93 xPrime[YY] += deltas[i % 4];
95 xPrime[ZZ] += deltas[i % 4];
98 std::fill(v_.begin(), v_.end(), RVec {0.0, 0.0, 0.0});
100 // Set up the topology.
101 const int settleType = 0;
102 mtop_.moltype.resize(1);
103 mtop_.molblock.resize(1);
104 mtop_.molblock[0].type = 0;
105 std::vector<int> &iatoms = mtop_.moltype[0].ilist[F_SETTLE].iatoms;
106 for (int i = 0; i < numSettles; ++i)
108 iatoms.push_back(settleType);
109 iatoms.push_back(i*atomsPerSettle_ + 0);
110 iatoms.push_back(i*atomsPerSettle_ + 1);
111 iatoms.push_back(i*atomsPerSettle_ + 2);
114 // Set up the SETTLE parameters.
116 iparams.settle.doh = dOH_;
117 iparams.settle.dhh = dHH_;
118 mtop_.ffparams.iparams.push_back(iparams);
120 // Set up the masses.
121 mtop_.moltype[0].atoms.atom = static_cast<t_atom*>(calloc(numSettles*atomsPerSettle_, sizeof(t_atom)));
122 mdatoms_.homenr = numSettles * atomsPerSettle_;
123 mdatoms_.massT = static_cast<real*>(calloc(mdatoms_.homenr, sizeof(real)));
124 mdatoms_.invmass = static_cast<real*>(calloc(mdatoms_.homenr, sizeof(real)));
125 for (int i = 0; i < numSettles; ++i)
127 mdatoms_.massT[i*atomsPerSettle_ + 0] = oxygenMass_;
128 mdatoms_.massT[i*atomsPerSettle_ + 1] = hydrogenMass_;
129 mdatoms_.massT[i*atomsPerSettle_ + 2] = hydrogenMass_;
131 mdatoms_.invmass[i*atomsPerSettle_ + 0] = 1.0/oxygenMass_;
132 mdatoms_.invmass[i*atomsPerSettle_ + 1] = 1.0/hydrogenMass_;
133 mdatoms_.invmass[i*atomsPerSettle_ + 2] = 1.0/hydrogenMass_;
135 mtop_.moltype[0].atoms.atom[i*atomsPerSettle_ + 0].m = oxygenMass_;
136 mtop_.moltype[0].atoms.atom[i*atomsPerSettle_ + 1].m = hydrogenMass_;
137 mtop_.moltype[0].atoms.atom[i*atomsPerSettle_ + 2].m = hydrogenMass_;
140 // Reshape some data so it can be directly used by the SETTLE constraints
141 ilist_ = { mtop_.moltype[0].ilist[F_SETTLE].size(), 0, mtop_.moltype[0].ilist[F_SETTLE].iatoms.data(), 0 };
142 idef_.il[F_SETTLE] = ilist_;
145 SettleTestData::~SettleTestData()
147 free(mdatoms_.massT);
148 free(mdatoms_.invmass);