2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
36 * \brief SHAKE and LINCS tests.
38 * \todo Better tests for virial are needed.
39 * \todo Tests for bigger systems to test threads synchronization,
40 * reduction, etc. on the GPU.
41 * \todo Tests for algorithms for derivatives.
42 * \todo Free-energy perturbation tests
44 * \author Artem Zhmurov <zhmurov@gmail.com>
45 * \ingroup module_mdlib
53 #include <unordered_map>
56 #include <gtest/gtest.h>
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/utility/stringutil.h"
63 #include "testutils/test_hardware_environment.h"
64 #include "testutils/testasserts.h"
66 #include "constrtestdata.h"
67 #include "constrtestrunners.h"
76 // Define the set of PBCs to run the test for
77 const std::vector
<t_pbc
> c_pbcs
= [] {
78 std::vector
<t_pbc
> pbcs
;
81 // Infinitely small box
82 matrix boxNone
= { { 0, 0, 0 }, { 0, 0, 0 }, { 0, 0, 0 } };
83 set_pbc(&pbc
, PbcType::No
, boxNone
);
84 pbcs
.emplace_back(pbc
);
87 matrix boxXyz
= { { 10.0, 0.0, 0.0 }, { 0.0, 20.0, 0.0 }, { 0.0, 0.0, 15.0 } };
88 set_pbc(&pbc
, PbcType::Xyz
, boxXyz
);
89 pbcs
.emplace_back(pbc
);
94 /*! \brief Test fixture for constraints.
96 * The fixture uses following test systems:
97 * 1. Two atoms, connected with one constraint (e.g. NH).
98 * 2. Three atoms, connected consequently with two constraints (e.g. CH2).
99 * 3. Three atoms, constrained to the fourth atom (e.g. CH3).
100 * 4. Four atoms, connected by two independent constraints.
101 * 5. Three atoms, connected by three constraints in a triangle
102 * (e.g. H2O with constrained H-O-H angle).
103 * 6. Four atoms, connected by three consequential constraints.
105 * For all systems, the final lengths of the constraints are tested against the
106 * reference values, the direction of each constraint is checked.
107 * Test also verifies that the center of mass has not been
108 * shifted by the constraints and that its velocity has not changed.
109 * For some systems, the value for scaled virial tensor is checked against
112 class ConstraintsTest
: public ::testing::TestWithParam
<t_pbc
>
116 * The test on the final length of constrained bonds.
118 * Goes through all the constraints and checks if the final length of all the constraints is
119 * equal to the target length with provided tolerance.
121 * \param[in] tolerance Allowed tolerance in final lengths.
122 * \param[in] testData Test data structure.
123 * \param[in] pbc Periodic boundary data.
125 static void checkConstrainsLength(FloatingPointTolerance tolerance
,
126 const ConstraintsTestData
& testData
,
130 // Test if all the constraints are satisfied
131 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
133 real r0
= testData
.constraintsR0_
.at(testData
.constraints_
.at(3 * c
));
134 int i
= testData
.constraints_
.at(3 * c
+ 1);
135 int j
= testData
.constraints_
.at(3 * c
+ 2);
138 if (pbc
.pbcType
== PbcType::Xyz
)
140 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
141 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
145 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
146 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
150 EXPECT_REAL_EQ_TOL(r0
, d1
, tolerance
) << gmx::formatString(
151 "rij = %f, which is not equal to r0 = %f for constraint #%zd, between atoms %d "
153 " (before constraining rij was %f).",
154 d1
, r0
, c
, i
, j
, d0
);
159 * The test on the final length of constrained bonds.
161 * Goes through all the constraints and checks if the direction of constraint has not changed
162 * by the algorithm (i.e. the constraints algorithm arrived to the solution that is closest
163 * to the initial system conformation).
165 * \param[in] testData Test data structure.
166 * \param[in] pbc Periodic boundary data.
168 static void checkConstrainsDirection(const ConstraintsTestData
& testData
, t_pbc pbc
)
171 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
173 int i
= testData
.constraints_
.at(3 * c
+ 1);
174 int j
= testData
.constraints_
.at(3 * c
+ 2);
176 if (pbc
.pbcType
== PbcType::Xyz
)
178 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
179 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
183 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
184 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
187 real dot
= xij0
.dot(xij1
);
189 EXPECT_GE(dot
, 0.0) << gmx::formatString(
190 "The constraint %zd changed direction. Constraining algorithm might have "
191 "returned the wrong root "
192 "of the constraints equation.",
198 * The test on the coordinates of the center of the mass (COM) of the system.
200 * Checks if the center of mass has not been shifted by the constraints. Note,
201 * that this test does not take into account the periodic boundary conditions.
202 * Hence it will not work should the constraints decide to move atoms across
205 * \param[in] tolerance Allowed tolerance in COM coordinates.
206 * \param[in] testData Test data structure.
208 static void checkCOMCoordinates(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
211 RVec
comPrime0({ 0.0, 0.0, 0.0 });
212 RVec
comPrime({ 0.0, 0.0, 0.0 });
213 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
215 comPrime0
+= testData
.masses_
[i
] * testData
.xPrime0_
[i
];
216 comPrime
+= testData
.masses_
[i
] * testData
.xPrime_
[i
];
219 comPrime0
/= testData
.numAtoms_
;
220 comPrime
/= testData
.numAtoms_
;
222 EXPECT_REAL_EQ_TOL(comPrime
[XX
], comPrime0
[XX
], tolerance
)
223 << "Center of mass was shifted by constraints in x-direction.";
224 EXPECT_REAL_EQ_TOL(comPrime
[YY
], comPrime0
[YY
], tolerance
)
225 << "Center of mass was shifted by constraints in y-direction.";
226 EXPECT_REAL_EQ_TOL(comPrime
[ZZ
], comPrime0
[ZZ
], tolerance
)
227 << "Center of mass was shifted by constraints in z-direction.";
231 * The test on the velocity of the center of the mass (COM) of the system.
233 * Checks if the velocity of the center of mass has not changed.
235 * \param[in] tolerance Allowed tolerance in COM velocity components.
236 * \param[in] testData Test data structure.
238 static void checkCOMVelocity(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
241 RVec
comV0({ 0.0, 0.0, 0.0 });
242 RVec
comV({ 0.0, 0.0, 0.0 });
243 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
245 comV0
+= testData
.masses_
[i
] * testData
.v0_
[i
];
246 comV
+= testData
.masses_
[i
] * testData
.v_
[i
];
248 comV0
/= testData
.numAtoms_
;
249 comV
/= testData
.numAtoms_
;
251 EXPECT_REAL_EQ_TOL(comV
[XX
], comV0
[XX
], tolerance
)
252 << "Velocity of the center of mass in x-direction has been changed by constraints.";
253 EXPECT_REAL_EQ_TOL(comV
[YY
], comV0
[YY
], tolerance
)
254 << "Velocity of the center of mass in y-direction has been changed by constraints.";
255 EXPECT_REAL_EQ_TOL(comV
[ZZ
], comV0
[ZZ
], tolerance
)
256 << "Velocity of the center of mass in z-direction has been changed by constraints.";
260 * The test of virial tensor.
262 * Checks if the values in the scaled virial tensor are equal to pre-computed values.
264 * \param[in] tolerance Tolerance for the tensor values.
265 * \param[in] testData Test data structure.
267 static void checkVirialTensor(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
269 for (int i
= 0; i
< DIM
; i
++)
271 for (int j
= 0; j
< DIM
; j
++)
273 EXPECT_REAL_EQ_TOL(testData
.virialScaledRef_
[i
][j
], testData
.virialScaled_
[i
][j
], tolerance
)
274 << gmx::formatString(
275 "Values in virial tensor at [%d][%d] are not within the "
276 "tolerance from reference value.",
281 //! Before any test is run, work out whether any compatible GPUs exist.
282 static std::vector
<std::unique_ptr
<IConstraintsTestRunner
>> getRunners()
284 std::vector
<std::unique_ptr
<IConstraintsTestRunner
>> runners
;
285 // Add runners for CPU versions of SHAKE and LINCS
286 runners
.emplace_back(std::make_unique
<ShakeConstraintsRunner
>());
287 runners
.emplace_back(std::make_unique
<LincsConstraintsRunner
>());
288 // If using CUDA, add runners for the GPU version of LINCS for each available GPU
291 for (const auto& testDevice
: getTestHardwareEnvironment()->getTestDeviceList())
293 runners
.emplace_back(std::make_unique
<LincsDeviceConstraintsRunner
>(*testDevice
));
300 TEST_P(ConstraintsTest
, SingleConstraint
)
302 std::string title
= "one constraint (e.g. OH)";
305 std::vector
<real
> masses
= { 1.0, 12.0 };
306 std::vector
<int> constraints
= { 0, 0, 1 };
307 std::vector
<real
> constraintsR0
= { 0.1 };
309 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
311 std::vector
<RVec
> x
= { { 0.0, oneTenthOverSqrtTwo
, 0.0 }, { oneTenthOverSqrtTwo
, 0.0, 0.0 } };
312 std::vector
<RVec
> xPrime
= { { 0.01, 0.08, 0.01 }, { 0.06, 0.01, -0.01 } };
313 std::vector
<RVec
> v
= { { 1.0, 2.0, 3.0 }, { 3.0, 2.0, 1.0 } };
315 tensor virialScaledRef
= { { -5.58e-04, 5.58e-04, 0.00e+00 },
316 { 5.58e-04, -5.58e-04, 0.00e+00 },
317 { 0.00e+00, 0.00e+00, 0.00e+00 } };
319 real shakeTolerance
= 0.0001;
320 gmx_bool shakeUseSOR
= false;
323 int lincslincsExpansionOrder
= 4;
324 real lincsWarnAngle
= 30.0;
326 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
327 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
328 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
329 lincslincsExpansionOrder
, lincsWarnAngle
);
331 t_pbc pbc
= GetParam();
333 // Cycle through all available runners
334 for (const auto& runner
: getRunners())
336 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
337 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
342 runner
->applyConstraints(testData
.get(), pbc
);
344 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
345 checkConstrainsDirection(*testData
, pbc
);
346 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
347 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
349 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
353 TEST_P(ConstraintsTest
, TwoDisjointConstraints
)
356 std::string title
= "two disjoint constraints";
358 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
359 std::vector
<int> constraints
= { 0, 0, 1, 1, 2, 3 };
360 std::vector
<real
> constraintsR0
= { 2.0, 1.0 };
363 std::vector
<RVec
> x
= {
364 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
367 std::vector
<RVec
> xPrime
= {
368 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
371 std::vector
<RVec
> v
= { { 0.0, 1.0, 0.0 }, { 1.0, 0.0, 0.0 }, { 0.0, 0.0, 1.0 }, { 0.0, 0.0, 0.0 } };
373 tensor virialScaledRef
= { { 3.3e-03, -1.7e-04, 5.6e-04 },
374 { -1.7e-04, 8.9e-06, -2.8e-05 },
375 { 5.6e-04, -2.8e-05, 8.9e-05 } };
377 real shakeTolerance
= 0.0001;
378 gmx_bool shakeUseSOR
= false;
381 int lincslincsExpansionOrder
= 4;
382 real lincsWarnAngle
= 30.0;
384 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
385 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
386 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
387 lincslincsExpansionOrder
, lincsWarnAngle
);
389 t_pbc pbc
= GetParam();
391 // Cycle through all available runners
392 for (const auto& runner
: getRunners())
394 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
395 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
400 runner
->applyConstraints(testData
.get(), pbc
);
402 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
403 checkConstrainsDirection(*testData
, pbc
);
404 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
405 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
407 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
411 TEST_P(ConstraintsTest
, ThreeSequentialConstraints
)
414 std::string title
= "three atoms, connected longitudinally (e.g. CH2)";
416 std::vector
<real
> masses
= { 1.0, 12.0, 16.0 };
417 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2 };
418 std::vector
<real
> constraintsR0
= { 0.1, 0.2 };
420 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
421 real twoTenthsOverSqrtThree
= 0.2_real
/ std::sqrt(3.0_real
);
423 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, oneTenthOverSqrtTwo
, 0.0 },
425 { twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
} };
427 std::vector
<RVec
> xPrime
= { { 0.08, 0.07, 0.01 }, { -0.02, 0.01, -0.02 }, { 0.10, 0.12, 0.11 } };
429 std::vector
<RVec
> v
= { { 1.0, 0.0, 0.0 }, { 0.0, 1.0, 0.0 }, { 0.0, 0.0, 1.0 } };
431 tensor virialScaledRef
= { { 4.14e-03, 4.14e-03, 3.31e-03 },
432 { 4.14e-03, 4.14e-03, 3.31e-03 },
433 { 3.31e-03, 3.31e-03, 3.31e-03 } };
435 real shakeTolerance
= 0.0001;
436 gmx_bool shakeUseSOR
= false;
439 int lincslincsExpansionOrder
= 4;
440 real lincsWarnAngle
= 30.0;
442 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
443 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
444 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
445 lincslincsExpansionOrder
, lincsWarnAngle
);
447 t_pbc pbc
= GetParam();
449 // Cycle through all available runners
450 for (const auto& runner
: getRunners())
452 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
453 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
458 runner
->applyConstraints(testData
.get(), pbc
);
460 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
461 checkConstrainsDirection(*testData
, pbc
);
462 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
463 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
465 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
469 TEST_P(ConstraintsTest
, ThreeConstraintsWithCentralAtom
)
472 std::string title
= "three atoms, connected to the central atom (e.g. CH3)";
474 std::vector
<real
> masses
= { 12.0, 1.0, 1.0, 1.0 };
475 std::vector
<int> constraints
= { 0, 0, 1, 0, 0, 2, 0, 0, 3 };
476 std::vector
<real
> constraintsR0
= { 0.1 };
479 std::vector
<RVec
> x
= {
480 { 0.00, 0.00, 0.00 }, { 0.10, 0.00, 0.00 }, { 0.00, -0.10, 0.00 }, { 0.00, 0.00, 0.10 }
483 std::vector
<RVec
> xPrime
= { { 0.004, 0.009, -0.010 },
484 { 0.110, -0.006, 0.003 },
485 { -0.007, -0.102, -0.007 },
486 { -0.005, 0.011, 0.102 } };
488 std::vector
<RVec
> v
= { { 1.0, 0.0, 0.0 }, { 1.0, 0.0, 0.0 }, { 1.0, 0.0, 0.0 }, { 1.0, 0.0, 0.0 } };
490 tensor virialScaledRef
= { { 7.14e-04, 0.00e+00, 0.00e+00 },
491 { 0.00e+00, 1.08e-03, 0.00e+00 },
492 { 0.00e+00, 0.00e+00, 1.15e-03 } };
494 real shakeTolerance
= 0.0001;
495 gmx_bool shakeUseSOR
= false;
498 int lincslincsExpansionOrder
= 4;
499 real lincsWarnAngle
= 30.0;
501 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
502 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
503 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
504 lincslincsExpansionOrder
, lincsWarnAngle
);
506 t_pbc pbc
= GetParam();
508 // Cycle through all available runners
509 for (const auto& runner
: getRunners())
511 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
512 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
517 runner
->applyConstraints(testData
.get(), pbc
);
519 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
520 checkConstrainsDirection(*testData
, pbc
);
521 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
522 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
524 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
528 TEST_P(ConstraintsTest
, FourSequentialConstraints
)
531 std::string title
= "four atoms, connected longitudinally";
533 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
534 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2, 2, 2, 3 };
535 std::vector
<real
> constraintsR0
= { 2.0, 1.0, 1.0 };
538 std::vector
<RVec
> x
= {
539 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
542 std::vector
<RVec
> xPrime
= {
543 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
546 std::vector
<RVec
> v
= { { 0.0, 0.0, 2.0 }, { 0.0, 0.0, 3.0 }, { 0.0, 0.0, -4.0 }, { 0.0, 0.0, -1.0 } };
548 tensor virialScaledRef
= { { 1.15e-01, -4.20e-03, 2.12e-02 },
549 { -4.20e-03, 1.70e-04, -6.41e-04 },
550 { 2.12e-02, -6.41e-04, 5.45e-03 } };
552 real shakeTolerance
= 0.0001;
553 gmx_bool shakeUseSOR
= false;
556 int lincslincsExpansionOrder
= 8;
557 real lincsWarnAngle
= 30.0;
559 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
560 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
561 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
562 lincslincsExpansionOrder
, lincsWarnAngle
);
564 t_pbc pbc
= GetParam();
566 // Cycle through all available runners
567 for (const auto& runner
: getRunners())
569 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
570 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
575 runner
->applyConstraints(testData
.get(), pbc
);
577 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
578 checkConstrainsDirection(*testData
, pbc
);
579 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
580 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
582 checkVirialTensor(absoluteTolerance(0.01), *testData
);
586 TEST_P(ConstraintsTest
, TriangleOfConstraints
)
589 std::string title
= "basic triangle (tree atoms, connected to each other)";
591 std::vector
<real
> masses
= { 1.0, 1.0, 1.0 };
592 std::vector
<int> constraints
= { 0, 0, 1, 2, 0, 2, 1, 1, 2 };
593 std::vector
<real
> constraintsR0
= { 0.1, 0.1, 0.1 };
595 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
597 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, 0.0, 0.0 },
598 { 0.0, oneTenthOverSqrtTwo
, 0.0 },
599 { 0.0, 0.0, oneTenthOverSqrtTwo
} };
601 std::vector
<RVec
> xPrime
= { { 0.09, -0.02, 0.01 }, { -0.02, 0.10, -0.02 }, { 0.03, -0.01, 0.07 } };
603 std::vector
<RVec
> v
= { { 1.0, 1.0, 1.0 }, { -2.0, -2.0, -2.0 }, { 1.0, 1.0, 1.0 } };
605 tensor virialScaledRef
= { { 6.00e-04, -1.61e-03, 1.01e-03 },
606 { -1.61e-03, 2.53e-03, -9.25e-04 },
607 { 1.01e-03, -9.25e-04, -8.05e-05 } };
609 real shakeTolerance
= 0.0001;
610 gmx_bool shakeUseSOR
= false;
613 int lincslincsExpansionOrder
= 4;
614 real lincsWarnAngle
= 30.0;
616 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
617 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
618 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
619 lincslincsExpansionOrder
, lincsWarnAngle
);
621 t_pbc pbc
= GetParam();
623 // Cycle through all available runners
624 for (const auto& runner
: getRunners())
626 SCOPED_TRACE(formatString("Testing %s with %s PBC using %s.", testData
->title_
.c_str(),
627 c_pbcTypeNames
[pbc
.pbcType
].c_str(), runner
->name().c_str()));
632 runner
->applyConstraints(testData
.get(), pbc
);
634 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
635 checkConstrainsDirection(*testData
, pbc
);
636 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
637 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
639 checkVirialTensor(absoluteTolerance(0.00001), *testData
);
643 INSTANTIATE_TEST_CASE_P(WithParameters
, ConstraintsTest
, ::testing::ValuesIn(c_pbcs
));