f717fc0207b7091c65f0e8a071d8687aed52e2b7
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / settle.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2018, 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 #include "gmxpre.h"
36
37 #include "gromacs/mdlib/settle.h"
38
39 #include <tuple>
40 #include <vector>
41
42 #include <gtest/gtest.h>
43
44 #include "gromacs/math/vec.h"
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/mdtypes/mdatom.h"
47 #include "gromacs/pbcutil/pbc.h"
48 #include "gromacs/topology/idef.h"
49 #include "gromacs/topology/ifunc.h"
50 #include "gromacs/topology/topology.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/utility/stringutil.h"
53 #include "gromacs/utility/unique_cptr.h"
54
55 #include "testutils/testasserts.h"
56
57 namespace gmx
58 {
59
60 namespace test
61 {
62
63 //! Database of 51 water atom input positions (DIM reals per atom, taken from spc216.gro) for use as test inputs.
64 const double g_positions[] = {
65     .130, -.041, -.291,
66     .120, -.056, -.192,
67     .044, -.005, -.327,
68     -.854, -.406, .477,
69     -.900, -.334, .425,
70     -.858, -.386, .575,
71     .351, -.061, .853,
72     .401, -.147, .859,
73     .416, .016, .850,
74     -.067, -.796, .873,
75     -.129, -.811, .797,
76     -.119, -.785, .958,
77     -.635, -.312, -.356,
78     -.629, -.389, -.292,
79     -.687, -.338, -.436,
80     .321, -.919, .242,
81     .403, -.880, .200,
82     .294, -1.001, .193,
83     -.404, .735, .728,
84     -.409, .670, .803,
85     -.324, .794, .741,
86     .461, -.596, -.135,
87     .411, -.595, -.221,
88     .398, -.614, -.059,
89     -.751, -.086, .237,
90     -.811, -.148, .287,
91     -.720, -.130, .152,
92     .202, .285, -.364,
93     .122, .345, -.377,
94     .192, .236, -.278,
95     -.230, -.485, .081,
96     -.262, -.391, .071,
97     -.306, -.548, .069,
98     .464, -.119, .323,
99     .497, -.080, .409,
100     .540, -.126, .258,
101     -.462, .107, .426,
102     -.486, .070, .336,
103     -.363, .123, .430,
104     .249, -.077, -.621,
105     .306, -.142, -.571,
106     .233, -.110, -.714,
107     -.922, -.164, .904,
108     -.842, -.221, .925,
109     -.971, -.204, .827,
110     .382, .700, .480,
111     .427, .610, .477,
112     .288, .689, .513,
113     .781, .264, -.113,
114     .848, .203, -.070,
115     .708, .283, -.048
116 };
117
118 //! Simple cubic simulation box to use in tests
119 matrix g_box = {{real(1.86206), 0, 0}, {0, real(1.86206), 0}, {0, 0, real(1.86206)}};
120
121 //! Convenience typedef
122 typedef std::tuple<int, bool, bool, bool> SettleTestParameters;
123
124 /*! \brief Test fixture for testing SETTLE position updates
125  *
126  * \todo This also tests that if the calling code requires velocities
127  * and virial updates, that those outputs do change, but does not test
128  * that those changes are correct.
129  *
130  * \todo Only no-PBC and cubic-PBC are tested here, but the correct
131  * function of the SIMD version of set_pbx_auic in all cases should be
132  * tested elsewhere.
133  */
134 class SettleTest : public ::testing::TestWithParam<SettleTestParameters>
135 {
136     public:
137         //! Updated water atom positions to constrain (DIM reals per atom)
138         std::vector<real> updatedPositions_;
139         //! Water atom velocities to constrain (DIM reals per atom)
140         std::vector<real> velocities_;
141         //! PBC option to test
142         t_pbc             pbcNone_;
143         //! PBC option to test
144         t_pbc             pbcXYZ_;
145
146         //! Constructor
147         SettleTest() :
148             updatedPositions_(std::begin(g_positions), std::end(g_positions)),
149             velocities_(updatedPositions_.size(), 0)
150         {
151             set_pbc(&pbcNone_, epbcNONE, g_box);
152             set_pbc(&pbcXYZ_, epbcXYZ, g_box);
153
154             // Perturb the atom positions, to appear like an
155             // "update," and where there is definitely constraining
156             // work to do.
157             for (size_t i = 0; i != updatedPositions_.size(); ++i)
158             {
159                 if (i % 4 == 0)
160                 {
161                     updatedPositions_[i] += 0.01;
162                 }
163                 else if (i % 4 == 1)
164                 {
165                     updatedPositions_[i] -= 0.01;
166                 }
167                 else if (i % 4 == 2)
168                 {
169                     updatedPositions_[i] += 0.02;
170                 }
171                 else if (i % 4 == 3)
172                 {
173                     updatedPositions_[i] -= 0.02;
174                 }
175             }
176         }
177 };
178
179 TEST_P(SettleTest, SatisfiesConstraints)
180 {
181     int  numSettles;
182     bool usePbc, useVelocities, calcVirial;
183     // Make some symbolic names for the parameter combination under
184     // test.
185     std::tie(numSettles, usePbc, useVelocities, calcVirial) = GetParam();
186
187     // Make a string that describes which parameter combination is
188     // being tested, to help make failing tests comprehensible.
189     std::string testDescription = formatString("while testing %d SETTLEs, %sPBC, %svelocities and %scalculating the virial",
190                                                numSettles,
191                                                usePbc ? "with " : "without ",
192                                                useVelocities ? "with " : "without ",
193                                                calcVirial ? "" : "not ");
194
195     const int settleType     = 0;
196     const int atomsPerSettle = NRAL(F_SETTLE);
197     ASSERT_LE(numSettles, updatedPositions_.size() / (atomsPerSettle * DIM)) << "cannot test that many SETTLEs " << testDescription;
198
199     // Set up the topology. We still have to make some raw pointers,
200     // but they are put into scope guards for automatic cleanup.
201     gmx_mtop_t                   *mtop;
202     snew(mtop, 1);
203     const unique_cptr<gmx_mtop_t> mtopGuard(mtop);
204     mtop->nmoltype = 1;
205     snew(mtop->moltype, mtop->nmoltype);
206     const unique_cptr<gmx_moltype_t> moltypeGuard(mtop->moltype);
207     mtop->nmolblock = 1;
208     snew(mtop->molblock, mtop->nmolblock);
209     const unique_cptr<gmx_molblock_t> molblockGuard(mtop->molblock);
210     mtop->molblock[0].type = 0;
211     std::vector<int>                  iatoms;
212     for (int i = 0; i < numSettles; ++i)
213     {
214         iatoms.push_back(settleType);
215         iatoms.push_back(i*atomsPerSettle+0);
216         iatoms.push_back(i*atomsPerSettle+1);
217         iatoms.push_back(i*atomsPerSettle+2);
218     }
219     mtop->moltype[0].ilist[F_SETTLE].iatoms = iatoms.data();
220     mtop->moltype[0].ilist[F_SETTLE].nr     = iatoms.size();
221
222     // Set up the SETTLE parameters.
223     mtop->ffparams.ntypes = 1;
224     snew(mtop->ffparams.iparams, mtop->ffparams.ntypes);
225     const unique_cptr<t_iparams> iparamsGuard(mtop->ffparams.iparams);
226     const real                   dOH = 0.09572;
227     const real                   dHH = 0.15139;
228     mtop->ffparams.iparams[settleType].settle.doh = dOH;
229     mtop->ffparams.iparams[settleType].settle.dhh = dHH;
230
231     // Set up the masses.
232     t_mdatoms         mdatoms;
233     std::vector<real> mass, massReciprocal;
234     const real        oxygenMass = 15.9994, hydrogenMass = 1.008;
235     for (int i = 0; i < numSettles; ++i)
236     {
237         mass.push_back(oxygenMass);
238         mass.push_back(hydrogenMass);
239         mass.push_back(hydrogenMass);
240         massReciprocal.push_back(1./oxygenMass);
241         massReciprocal.push_back(1./hydrogenMass);
242         massReciprocal.push_back(1./hydrogenMass);
243     }
244     mdatoms.massT   = mass.data();
245     mdatoms.invmass = massReciprocal.data();
246     mdatoms.homenr  = numSettles * atomsPerSettle;
247
248     // Finally make the settle data structures
249     gmx_settledata_t settled = settle_init(mtop);
250     settle_set_constraints(settled, &mtop->moltype[0].ilist[F_SETTLE], &mdatoms);
251
252     // Copy the original positions from the array of doubles to a vector of reals
253     std::vector<real> startingPositions(std::begin(g_positions), std::end(g_positions));
254
255     // Run the test
256     bool       errorOccured;
257     tensor     virial             = {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}};
258     int        numThreads         = 1, threadIndex = 0;
259     const real reciprocalTimeStep = 1.0/0.002;
260     csettle(settled, numThreads, threadIndex,
261             usePbc ? &pbcXYZ_ : &pbcNone_,
262             startingPositions.data(), updatedPositions_.data(), reciprocalTimeStep,
263             useVelocities ? velocities_.data() : nullptr,
264             calcVirial, virial, &errorOccured);
265     settle_free(settled);
266     EXPECT_FALSE(errorOccured) << testDescription;
267
268     // The necessary tolerances for the test to pass were determined
269     // empirically. This isn't nice, but the required behaviour that
270     // SETTLE produces constrained coordinates consistent with
271     // sensible sampling needs to be tested at a much higher level.
272     FloatingPointTolerance tolerance =
273         relativeToleranceAsPrecisionDependentUlp(dOH*dOH, 80, 380);
274
275     // Verify the updated coordinates match the requirements
276     int positionIndex = 0, velocityIndex = 0;
277     for (int i = 0; i < numSettles; ++i)
278     {
279         rvec positionO  {
280             updatedPositions_[positionIndex++], updatedPositions_[positionIndex++], updatedPositions_[positionIndex++]
281         };
282         rvec positionH1 {
283             updatedPositions_[positionIndex++], updatedPositions_[positionIndex++], updatedPositions_[positionIndex++]
284         };
285         rvec positionH2 {
286             updatedPositions_[positionIndex++], updatedPositions_[positionIndex++], updatedPositions_[positionIndex++]
287         };
288
289         EXPECT_REAL_EQ_TOL(dOH*dOH, distance2(positionO, positionH1), tolerance) << formatString("for water %d ", i) << testDescription;
290         EXPECT_REAL_EQ_TOL(dOH*dOH, distance2(positionO, positionH2), tolerance) << formatString("for water %d ", i) << testDescription;
291         EXPECT_REAL_EQ_TOL(dHH*dHH, distance2(positionH1, positionH2), tolerance) << formatString("for water %d ", i) << testDescription;
292
293         // This merely tests whether the velocities were
294         // updated from the starting values of zero (or not),
295         // but not whether the update was correct.
296         for (int j = 0; j < atomsPerSettle * DIM; ++j, ++velocityIndex)
297         {
298             EXPECT_TRUE(useVelocities == (0. != velocities_[velocityIndex])) << formatString("for water %d velocity coordinate %d ", i, j) << testDescription;
299         }
300     }
301
302     // This merely tests whether the viral was updated from
303     // the starting values of zero (or not), but not whether
304     // any update was correct.
305     for (int d = 0; d < DIM; ++d)
306     {
307         for (int dd = 0; dd < DIM; ++dd)
308         {
309             EXPECT_TRUE(calcVirial == (0. != virial[d][dd])) << formatString("for virial component[%d][%d] ", d, dd) << testDescription;
310         }
311     }
312 }
313
314 // Scan the full Cartesian product of numbers of SETTLE interactions
315 // (4 and 17 are chosen to test cases that do and do not match
316 // hardware SIMD widths), and whether or not we use PBC, velocities or
317 // calculate the virial contribution.
318 INSTANTIATE_TEST_CASE_P(WithParameters, SettleTest,
319                             ::testing::Combine(::testing::Values(1, 4, 7),
320                                                    ::testing::Bool(),
321                                                    ::testing::Bool(),
322                                                    ::testing::Bool()));
323
324 } // namespace
325 } // namespace