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 /*! \brief Test fixture for constraints.
78 * The fixture uses following test systems:
79 * 1. Two atoms, connected with one constraint (e.g. NH).
80 * 2. Three atoms, connected consequently with two constraints (e.g. CH2).
81 * 3. Three atoms, constrained to the fourth atom (e.g. CH3).
82 * 4. Four atoms, connected by two independent constraints.
83 * 5. Three atoms, connected by three constraints in a triangle
84 * (e.g. H2O with constrained H-O-H angle).
85 * 6. Four atoms, connected by three consequential constraints.
87 * For all systems, the final lengths of the constraints are tested against the
88 * reference values, the direction of each constraint is checked.
89 * Test also verifies that the center of mass has not been
90 * shifted by the constraints and that its velocity has not changed.
91 * For some systems, the value for scaled virial tensor is checked against
94 class ConstraintsTest
: public ::testing::TestWithParam
<std::string
>
98 std::unordered_map
<std::string
, t_pbc
> pbcs_
;
100 /*! \brief Test setup function.
102 * Setting up the pbcs and algorithms. Note, that corresponding string keywords
103 * have to be explicitly added at the end of this file when the tests are called.
106 void SetUp() override
110 // PBC initialization
114 // Infinitely small box
115 matrix boxNone
= { { 0, 0, 0 }, { 0, 0, 0 }, { 0, 0, 0 } };
116 set_pbc(&pbc
, PbcType::No
, boxNone
);
117 pbcs_
["PBCNone"] = pbc
;
120 matrix boxXyz
= { { 10.0, 0.0, 0.0 }, { 0.0, 20.0, 0.0 }, { 0.0, 0.0, 15.0 } };
121 set_pbc(&pbc
, PbcType::Xyz
, boxXyz
);
122 pbcs_
["PBCXYZ"] = pbc
;
126 * The test on the final length of constrained bonds.
128 * Goes through all the constraints and checks if the final length of all the constraints is
129 * equal to the target length with provided tolerance.
131 * \param[in] tolerance Allowed tolerance in final lengths.
132 * \param[in] testData Test data structure.
133 * \param[in] pbc Periodic boundary data.
135 static void checkConstrainsLength(FloatingPointTolerance tolerance
,
136 const ConstraintsTestData
& testData
,
140 // Test if all the constraints are satisfied
141 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
143 real r0
= testData
.constraintsR0_
.at(testData
.constraints_
.at(3 * c
));
144 int i
= testData
.constraints_
.at(3 * c
+ 1);
145 int j
= testData
.constraints_
.at(3 * c
+ 2);
148 if (pbc
.pbcType
== PbcType::Xyz
)
150 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
151 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
155 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
156 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
160 EXPECT_REAL_EQ_TOL(r0
, d1
, tolerance
) << gmx::formatString(
161 "rij = %f, which is not equal to r0 = %f for constraint #%zd, between atoms %d "
163 " (before constraining rij was %f).",
164 d1
, r0
, c
, i
, j
, d0
);
169 * The test on the final length of constrained bonds.
171 * Goes through all the constraints and checks if the direction of constraint has not changed
172 * by the algorithm (i.e. the constraints algorithm arrived to the solution that is closest
173 * to the initial system conformation).
175 * \param[in] testData Test data structure.
176 * \param[in] pbc Periodic boundary data.
178 static void checkConstrainsDirection(const ConstraintsTestData
& testData
, t_pbc pbc
)
181 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
183 int i
= testData
.constraints_
.at(3 * c
+ 1);
184 int j
= testData
.constraints_
.at(3 * c
+ 2);
186 if (pbc
.pbcType
== PbcType::Xyz
)
188 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
189 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
193 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
194 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
197 real dot
= xij0
.dot(xij1
);
199 EXPECT_GE(dot
, 0.0) << gmx::formatString(
200 "The constraint %zd changed direction. Constraining algorithm might have "
201 "returned the wrong root "
202 "of the constraints equation.",
208 * The test on the coordinates of the center of the mass (COM) of the system.
210 * Checks if the center of mass has not been shifted by the constraints. Note,
211 * that this test does not take into account the periodic boundary conditions.
212 * Hence it will not work should the constraints decide to move atoms across
215 * \param[in] tolerance Allowed tolerance in COM coordinates.
216 * \param[in] testData Test data structure.
218 static void checkCOMCoordinates(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
221 RVec
comPrime0({ 0.0, 0.0, 0.0 });
222 RVec
comPrime({ 0.0, 0.0, 0.0 });
223 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
225 comPrime0
+= testData
.masses_
[i
] * testData
.xPrime0_
[i
];
226 comPrime
+= testData
.masses_
[i
] * testData
.xPrime_
[i
];
229 comPrime0
/= testData
.numAtoms_
;
230 comPrime
/= testData
.numAtoms_
;
232 EXPECT_REAL_EQ_TOL(comPrime
[XX
], comPrime0
[XX
], tolerance
)
233 << "Center of mass was shifted by constraints in x-direction.";
234 EXPECT_REAL_EQ_TOL(comPrime
[YY
], comPrime0
[YY
], tolerance
)
235 << "Center of mass was shifted by constraints in y-direction.";
236 EXPECT_REAL_EQ_TOL(comPrime
[ZZ
], comPrime0
[ZZ
], tolerance
)
237 << "Center of mass was shifted by constraints in z-direction.";
241 * The test on the velocity of the center of the mass (COM) of the system.
243 * Checks if the velocity of the center of mass has not changed.
245 * \param[in] tolerance Allowed tolerance in COM velocity components.
246 * \param[in] testData Test data structure.
248 static void checkCOMVelocity(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
251 RVec
comV0({ 0.0, 0.0, 0.0 });
252 RVec
comV({ 0.0, 0.0, 0.0 });
253 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
255 comV0
+= testData
.masses_
[i
] * testData
.v0_
[i
];
256 comV
+= testData
.masses_
[i
] * testData
.v_
[i
];
258 comV0
/= testData
.numAtoms_
;
259 comV
/= testData
.numAtoms_
;
261 EXPECT_REAL_EQ_TOL(comV
[XX
], comV0
[XX
], tolerance
)
262 << "Velocity of the center of mass in x-direction has been changed by constraints.";
263 EXPECT_REAL_EQ_TOL(comV
[YY
], comV0
[YY
], tolerance
)
264 << "Velocity of the center of mass in y-direction has been changed by constraints.";
265 EXPECT_REAL_EQ_TOL(comV
[ZZ
], comV0
[ZZ
], tolerance
)
266 << "Velocity of the center of mass in z-direction has been changed by constraints.";
270 * The test of virial tensor.
272 * Checks if the values in the scaled virial tensor are equal to pre-computed values.
274 * \param[in] tolerance Tolerance for the tensor values.
275 * \param[in] testData Test data structure.
277 static void checkVirialTensor(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
279 for (int i
= 0; i
< DIM
; i
++)
281 for (int j
= 0; j
< DIM
; j
++)
283 EXPECT_REAL_EQ_TOL(testData
.virialScaledRef_
[i
][j
], testData
.virialScaled_
[i
][j
], tolerance
)
284 << gmx::formatString(
285 "Values in virial tensor at [%d][%d] are not within the "
286 "tolerance from reference value.",
291 //! Before any test is run, work out whether any compatible GPUs exist.
292 static std::vector
<std::unique_ptr
<IConstraintsTestRunner
>> getRunners()
294 std::vector
<std::unique_ptr
<IConstraintsTestRunner
>> runners
;
295 // Add runners for CPU versions of SHAKE and LINCS
296 runners
.emplace_back(std::make_unique
<ShakeConstraintsRunner
>());
297 runners
.emplace_back(std::make_unique
<LincsConstraintsRunner
>());
298 // If using CUDA, add runners for the GPU version of LINCS for each available GPU
301 for (const auto& testDevice
: getTestHardwareEnvironment()->getTestDeviceList())
303 runners
.emplace_back(std::make_unique
<LincsDeviceConstraintsRunner
>(*testDevice
));
310 TEST_P(ConstraintsTest
, SingleConstraint
)
312 std::string title
= "one constraint (e.g. OH)";
315 std::vector
<real
> masses
= { 1.0, 12.0 };
316 std::vector
<int> constraints
= { 0, 0, 1 };
317 std::vector
<real
> constraintsR0
= { 0.1 };
319 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
321 std::vector
<RVec
> x
= { { 0.0, oneTenthOverSqrtTwo
, 0.0 }, { oneTenthOverSqrtTwo
, 0.0, 0.0 } };
322 std::vector
<RVec
> xPrime
= { { 0.01, 0.08, 0.01 }, { 0.06, 0.01, -0.01 } };
323 std::vector
<RVec
> v
= { { 1.0, 2.0, 3.0 }, { 3.0, 2.0, 1.0 } };
325 tensor virialScaledRef
= { { -5.58e-04, 5.58e-04, 0.00e+00 },
326 { 5.58e-04, -5.58e-04, 0.00e+00 },
327 { 0.00e+00, 0.00e+00, 0.00e+00 } };
329 real shakeTolerance
= 0.0001;
330 gmx_bool shakeUseSOR
= false;
333 int lincslincsExpansionOrder
= 4;
334 real lincsWarnAngle
= 30.0;
336 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
337 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
338 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
339 lincslincsExpansionOrder
, lincsWarnAngle
);
341 std::string pbcName
= GetParam();
342 t_pbc pbc
= pbcs_
.at(pbcName
);
344 // Cycle through all available runners
345 for (const auto& runner
: getRunners())
347 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
348 pbcName
.c_str(), runner
->name().c_str()));
353 runner
->applyConstraints(testData
.get(), pbc
);
355 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
356 checkConstrainsDirection(*testData
, pbc
);
357 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
358 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
360 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
364 TEST_P(ConstraintsTest
, TwoDisjointConstraints
)
367 std::string title
= "two disjoint constraints";
369 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
370 std::vector
<int> constraints
= { 0, 0, 1, 1, 2, 3 };
371 std::vector
<real
> constraintsR0
= { 2.0, 1.0 };
374 std::vector
<RVec
> x
= {
375 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
378 std::vector
<RVec
> xPrime
= {
379 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
382 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 } };
384 tensor virialScaledRef
= { { 3.3e-03, -1.7e-04, 5.6e-04 },
385 { -1.7e-04, 8.9e-06, -2.8e-05 },
386 { 5.6e-04, -2.8e-05, 8.9e-05 } };
388 real shakeTolerance
= 0.0001;
389 gmx_bool shakeUseSOR
= false;
392 int lincslincsExpansionOrder
= 4;
393 real lincsWarnAngle
= 30.0;
395 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
396 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
397 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
398 lincslincsExpansionOrder
, lincsWarnAngle
);
400 std::string pbcName
= GetParam();
401 t_pbc pbc
= pbcs_
.at(pbcName
);
403 // Cycle through all available runners
404 for (const auto& runner
: getRunners())
406 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
407 pbcName
.c_str(), runner
->name().c_str()));
412 runner
->applyConstraints(testData
.get(), pbc
);
414 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
415 checkConstrainsDirection(*testData
, pbc
);
416 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
417 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
419 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
423 TEST_P(ConstraintsTest
, ThreeSequentialConstraints
)
426 std::string title
= "three atoms, connected longitudinally (e.g. CH2)";
428 std::vector
<real
> masses
= { 1.0, 12.0, 16.0 };
429 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2 };
430 std::vector
<real
> constraintsR0
= { 0.1, 0.2 };
432 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
433 real twoTenthsOverSqrtThree
= 0.2_real
/ std::sqrt(3.0_real
);
435 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, oneTenthOverSqrtTwo
, 0.0 },
437 { twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
} };
439 std::vector
<RVec
> xPrime
= { { 0.08, 0.07, 0.01 }, { -0.02, 0.01, -0.02 }, { 0.10, 0.12, 0.11 } };
441 std::vector
<RVec
> v
= { { 1.0, 0.0, 0.0 }, { 0.0, 1.0, 0.0 }, { 0.0, 0.0, 1.0 } };
443 tensor virialScaledRef
= { { 4.14e-03, 4.14e-03, 3.31e-03 },
444 { 4.14e-03, 4.14e-03, 3.31e-03 },
445 { 3.31e-03, 3.31e-03, 3.31e-03 } };
447 real shakeTolerance
= 0.0001;
448 gmx_bool shakeUseSOR
= false;
451 int lincslincsExpansionOrder
= 4;
452 real lincsWarnAngle
= 30.0;
454 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
455 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
456 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
457 lincslincsExpansionOrder
, lincsWarnAngle
);
459 std::string pbcName
= GetParam();
460 t_pbc pbc
= pbcs_
.at(pbcName
);
462 // Cycle through all available runners
463 for (const auto& runner
: getRunners())
465 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
466 pbcName
.c_str(), runner
->name().c_str()));
471 runner
->applyConstraints(testData
.get(), pbc
);
473 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
474 checkConstrainsDirection(*testData
, pbc
);
475 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
476 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
478 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
482 TEST_P(ConstraintsTest
, ThreeConstraintsWithCentralAtom
)
485 std::string title
= "three atoms, connected to the central atom (e.g. CH3)";
487 std::vector
<real
> masses
= { 12.0, 1.0, 1.0, 1.0 };
488 std::vector
<int> constraints
= { 0, 0, 1, 0, 0, 2, 0, 0, 3 };
489 std::vector
<real
> constraintsR0
= { 0.1 };
492 std::vector
<RVec
> x
= {
493 { 0.00, 0.00, 0.00 }, { 0.10, 0.00, 0.00 }, { 0.00, -0.10, 0.00 }, { 0.00, 0.00, 0.10 }
496 std::vector
<RVec
> xPrime
= { { 0.004, 0.009, -0.010 },
497 { 0.110, -0.006, 0.003 },
498 { -0.007, -0.102, -0.007 },
499 { -0.005, 0.011, 0.102 } };
501 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 } };
503 tensor virialScaledRef
= { { 7.14e-04, 0.00e+00, 0.00e+00 },
504 { 0.00e+00, 1.08e-03, 0.00e+00 },
505 { 0.00e+00, 0.00e+00, 1.15e-03 } };
507 real shakeTolerance
= 0.0001;
508 gmx_bool shakeUseSOR
= false;
511 int lincslincsExpansionOrder
= 4;
512 real lincsWarnAngle
= 30.0;
514 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
515 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
516 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
517 lincslincsExpansionOrder
, lincsWarnAngle
);
519 std::string pbcName
= GetParam();
520 t_pbc pbc
= pbcs_
.at(pbcName
);
522 // Cycle through all available runners
523 for (const auto& runner
: getRunners())
525 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
526 pbcName
.c_str(), runner
->name().c_str()));
531 runner
->applyConstraints(testData
.get(), pbc
);
533 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
534 checkConstrainsDirection(*testData
, pbc
);
535 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
536 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
538 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
542 TEST_P(ConstraintsTest
, FourSequentialConstraints
)
545 std::string title
= "four atoms, connected longitudinally";
547 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
548 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2, 2, 2, 3 };
549 std::vector
<real
> constraintsR0
= { 2.0, 1.0, 1.0 };
552 std::vector
<RVec
> x
= {
553 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
556 std::vector
<RVec
> xPrime
= {
557 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
560 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 } };
562 tensor virialScaledRef
= { { 1.15e-01, -4.20e-03, 2.12e-02 },
563 { -4.20e-03, 1.70e-04, -6.41e-04 },
564 { 2.12e-02, -6.41e-04, 5.45e-03 } };
566 real shakeTolerance
= 0.0001;
567 gmx_bool shakeUseSOR
= false;
570 int lincslincsExpansionOrder
= 8;
571 real lincsWarnAngle
= 30.0;
573 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
574 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
575 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
576 lincslincsExpansionOrder
, lincsWarnAngle
);
578 std::string pbcName
= GetParam();
579 t_pbc pbc
= pbcs_
.at(pbcName
);
581 // Cycle through all available runners
582 for (const auto& runner
: getRunners())
584 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
585 pbcName
.c_str(), runner
->name().c_str()));
590 runner
->applyConstraints(testData
.get(), pbc
);
592 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
593 checkConstrainsDirection(*testData
, pbc
);
594 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
595 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
597 checkVirialTensor(absoluteTolerance(0.01), *testData
);
601 TEST_P(ConstraintsTest
, TriangleOfConstraints
)
604 std::string title
= "basic triangle (tree atoms, connected to each other)";
606 std::vector
<real
> masses
= { 1.0, 1.0, 1.0 };
607 std::vector
<int> constraints
= { 0, 0, 1, 2, 0, 2, 1, 1, 2 };
608 std::vector
<real
> constraintsR0
= { 0.1, 0.1, 0.1 };
610 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
612 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, 0.0, 0.0 },
613 { 0.0, oneTenthOverSqrtTwo
, 0.0 },
614 { 0.0, 0.0, oneTenthOverSqrtTwo
} };
616 std::vector
<RVec
> xPrime
= { { 0.09, -0.02, 0.01 }, { -0.02, 0.10, -0.02 }, { 0.03, -0.01, 0.07 } };
618 std::vector
<RVec
> v
= { { 1.0, 1.0, 1.0 }, { -2.0, -2.0, -2.0 }, { 1.0, 1.0, 1.0 } };
620 tensor virialScaledRef
= { { 6.00e-04, -1.61e-03, 1.01e-03 },
621 { -1.61e-03, 2.53e-03, -9.25e-04 },
622 { 1.01e-03, -9.25e-04, -8.05e-05 } };
624 real shakeTolerance
= 0.0001;
625 gmx_bool shakeUseSOR
= false;
628 int lincslincsExpansionOrder
= 4;
629 real lincsWarnAngle
= 30.0;
631 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
632 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
633 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
634 lincslincsExpansionOrder
, lincsWarnAngle
);
636 std::string pbcName
= GetParam();
637 t_pbc pbc
= pbcs_
.at(pbcName
);
639 // Cycle through all available runners
640 for (const auto& runner
: getRunners())
642 SCOPED_TRACE(formatString("Testing %s with %s using %s.", testData
->title_
.c_str(),
643 pbcName
.c_str(), runner
->name().c_str()));
648 runner
->applyConstraints(testData
.get(), pbc
);
650 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
651 checkConstrainsDirection(*testData
, pbc
);
652 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
653 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
655 checkVirialTensor(absoluteTolerance(0.00001), *testData
);
660 INSTANTIATE_TEST_CASE_P(WithParameters
, ConstraintsTest
, ::testing::Values("PBCNone", "PBCXYZ"));