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/testasserts.h"
65 #include "constrtestdata.h"
66 #include "constrtestrunners.h"
75 /*! \brief The two-dimensional parameter space for test.
77 * The test will run for all possible combinations of accessible
79 * 1. PBC setup ("PBCNONE" or "PBCXYZ")
80 * 2. The algorithm ("SHAKE", "LINCS" or "LINCS_GPU").
82 typedef std::tuple
<std::string
, std::string
> ConstraintsTestParameters
;
84 //! Names of all availible runners
85 std::vector
<std::string
> runnersNames
;
87 //! Method that fills and returns runnersNames to the test macros.
88 std::vector
<std::string
> getRunnersNames()
90 runnersNames
.emplace_back("SHAKE");
91 runnersNames
.emplace_back("LINCS");
92 if (GMX_GPU
== GMX_GPU_CUDA
&& canComputeOnGpu())
94 runnersNames
.emplace_back("LINCS_CUDA");
99 /*! \brief Test fixture for constraints.
101 * The fixture uses following test systems:
102 * 1. Two atoms, connected with one constraint (e.g. NH).
103 * 2. Three atoms, connected consequently with two constraints (e.g. CH2).
104 * 3. Three atoms, constrained to the fourth atom (e.g. CH3).
105 * 4. Four atoms, connected by two independent constraints.
106 * 5. Three atoms, connected by three constraints in a triangle
107 * (e.g. H2O with constrained H-O-H angle).
108 * 6. Four atoms, connected by three consequential constraints.
110 * For all systems, the final lengths of the constraints are tested against the
111 * reference values, the direction of each constraint is checked.
112 * Test also verifies that the center of mass has not been
113 * shifted by the constraints and that its velocity has not changed.
114 * For some systems, the value for scaled virial tensor is checked against
117 class ConstraintsTest
: public ::testing::TestWithParam
<ConstraintsTestParameters
>
121 std::unordered_map
<std::string
, t_pbc
> pbcs_
;
122 //! Algorithms (SHAKE and LINCS)
123 std::unordered_map
<std::string
, void (*)(ConstraintsTestData
* testData
, t_pbc pbc
)> algorithms_
;
125 /*! \brief Test setup function.
127 * Setting up the pbcs and algorithms. Note, that corresponding string keywords
128 * have to be explicitly added at the end of this file when the tests are called.
131 void SetUp() override
135 // PBC initialization
139 // Infinitely small box
140 matrix boxNone
= { { 0, 0, 0 }, { 0, 0, 0 }, { 0, 0, 0 } };
141 set_pbc(&pbc
, PbcType::No
, boxNone
);
142 pbcs_
["PBCNone"] = pbc
;
145 matrix boxXyz
= { { 10.0, 0.0, 0.0 }, { 0.0, 20.0, 0.0 }, { 0.0, 0.0, 15.0 } };
146 set_pbc(&pbc
, PbcType::Xyz
, boxXyz
);
147 pbcs_
["PBCXYZ"] = pbc
;
153 algorithms_
["SHAKE"] = applyShake
;
155 algorithms_
["LINCS"] = applyLincs
;
156 // LINCS using CUDA (will only be called if CUDA is available)
157 algorithms_
["LINCS_CUDA"] = applyLincsCuda
;
161 * The test on the final length of constrained bonds.
163 * Goes through all the constraints and checks if the final length of all the constraints is
164 * equal to the target length with provided tolerance.
166 * \param[in] tolerance Allowed tolerance in final lengths.
167 * \param[in] testData Test data structure.
168 * \param[in] pbc Periodic boundary data.
170 void checkConstrainsLength(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
, t_pbc pbc
)
173 // Test if all the constraints are satisfied
174 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
176 real r0
= testData
.constraintsR0_
.at(testData
.constraints_
.at(3 * c
));
177 int i
= testData
.constraints_
.at(3 * c
+ 1);
178 int j
= testData
.constraints_
.at(3 * c
+ 2);
181 if (pbc
.pbcType
== PbcType::Xyz
)
183 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
184 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
188 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
189 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
193 EXPECT_REAL_EQ_TOL(r0
, d1
, tolerance
) << gmx::formatString(
194 "rij = %f, which is not equal to r0 = %f for constraint #%zd, between atoms %d "
196 " (before constraining rij was %f).",
197 d1
, r0
, c
, i
, j
, d0
);
202 * The test on the final length of constrained bonds.
204 * Goes through all the constraints and checks if the direction of constraint has not changed
205 * by the algorithm (i.e. the constraints algorithm arrived to the solution that is closest
206 * to the initial system conformation).
208 * \param[in] testData Test data structure.
209 * \param[in] pbc Periodic boundary data.
211 void checkConstrainsDirection(const ConstraintsTestData
& testData
, t_pbc pbc
)
214 for (index c
= 0; c
< ssize(testData
.constraints_
) / 3; c
++)
216 int i
= testData
.constraints_
.at(3 * c
+ 1);
217 int j
= testData
.constraints_
.at(3 * c
+ 2);
219 if (pbc
.pbcType
== PbcType::Xyz
)
221 pbc_dx_aiuc(&pbc
, testData
.x_
[i
], testData
.x_
[j
], xij0
);
222 pbc_dx_aiuc(&pbc
, testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
226 rvec_sub(testData
.x_
[i
], testData
.x_
[j
], xij0
);
227 rvec_sub(testData
.xPrime_
[i
], testData
.xPrime_
[j
], xij1
);
230 real dot
= xij0
.dot(xij1
);
232 EXPECT_GE(dot
, 0.0) << gmx::formatString(
233 "The constraint %zd changed direction. Constraining algorithm might have "
234 "returned the wrong root "
235 "of the constraints equation.",
241 * The test on the coordinates of the center of the mass (COM) of the system.
243 * Checks if the center of mass has not been shifted by the constraints. Note,
244 * that this test does not take into account the periodic boundary conditions.
245 * Hence it will not work should the constraints decide to move atoms across
248 * \param[in] tolerance Allowed tolerance in COM coordinates.
249 * \param[in] testData Test data structure.
251 void checkCOMCoordinates(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
254 RVec
comPrime0({ 0.0, 0.0, 0.0 });
255 RVec
comPrime({ 0.0, 0.0, 0.0 });
256 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
258 comPrime0
+= testData
.masses_
[i
] * testData
.xPrime0_
[i
];
259 comPrime
+= testData
.masses_
[i
] * testData
.xPrime_
[i
];
262 comPrime0
/= testData
.numAtoms_
;
263 comPrime
/= testData
.numAtoms_
;
265 EXPECT_REAL_EQ_TOL(comPrime
[XX
], comPrime0
[XX
], tolerance
)
266 << "Center of mass was shifted by constraints in x-direction.";
267 EXPECT_REAL_EQ_TOL(comPrime
[YY
], comPrime0
[YY
], tolerance
)
268 << "Center of mass was shifted by constraints in y-direction.";
269 EXPECT_REAL_EQ_TOL(comPrime
[ZZ
], comPrime0
[ZZ
], tolerance
)
270 << "Center of mass was shifted by constraints in z-direction.";
274 * The test on the velocity of the center of the mass (COM) of the system.
276 * Checks if the velocity of the center of mass has not changed.
278 * \param[in] tolerance Allowed tolerance in COM velocity components.
279 * \param[in] testData Test data structure.
281 void checkCOMVelocity(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
284 RVec
comV0({ 0.0, 0.0, 0.0 });
285 RVec
comV({ 0.0, 0.0, 0.0 });
286 for (int i
= 0; i
< testData
.numAtoms_
; i
++)
288 comV0
+= testData
.masses_
[i
] * testData
.v0_
[i
];
289 comV
+= testData
.masses_
[i
] * testData
.v_
[i
];
291 comV0
/= testData
.numAtoms_
;
292 comV
/= testData
.numAtoms_
;
294 EXPECT_REAL_EQ_TOL(comV
[XX
], comV0
[XX
], tolerance
)
295 << "Velocity of the center of mass in x-direction has been changed by constraints.";
296 EXPECT_REAL_EQ_TOL(comV
[YY
], comV0
[YY
], tolerance
)
297 << "Velocity of the center of mass in y-direction has been changed by constraints.";
298 EXPECT_REAL_EQ_TOL(comV
[ZZ
], comV0
[ZZ
], tolerance
)
299 << "Velocity of the center of mass in z-direction has been changed by constraints.";
303 * The test of virial tensor.
305 * Checks if the values in the scaled virial tensor are equal to pre-computed values.
307 * \param[in] tolerance Tolerance for the tensor values.
308 * \param[in] testData Test data structure.
310 void checkVirialTensor(FloatingPointTolerance tolerance
, const ConstraintsTestData
& testData
)
312 for (int i
= 0; i
< DIM
; i
++)
314 for (int j
= 0; j
< DIM
; j
++)
316 EXPECT_REAL_EQ_TOL(testData
.virialScaledRef_
[i
][j
], testData
.virialScaled_
[i
][j
], tolerance
)
317 << gmx::formatString(
318 "Values in virial tensor at [%d][%d] are not within the "
319 "tolerance from reference value.",
326 TEST_P(ConstraintsTest
, SingleConstraint
)
328 std::string title
= "one constraint (e.g. OH)";
331 std::vector
<real
> masses
= { 1.0, 12.0 };
332 std::vector
<int> constraints
= { 0, 0, 1 };
333 std::vector
<real
> constraintsR0
= { 0.1 };
335 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
337 std::vector
<RVec
> x
= { { 0.0, oneTenthOverSqrtTwo
, 0.0 }, { oneTenthOverSqrtTwo
, 0.0, 0.0 } };
338 std::vector
<RVec
> xPrime
= { { 0.01, 0.08, 0.01 }, { 0.06, 0.01, -0.01 } };
339 std::vector
<RVec
> v
= { { 1.0, 2.0, 3.0 }, { 3.0, 2.0, 1.0 } };
341 tensor virialScaledRef
= { { -5.58e-04, 5.58e-04, 0.00e+00 },
342 { 5.58e-04, -5.58e-04, 0.00e+00 },
343 { 0.00e+00, 0.00e+00, 0.00e+00 } };
345 real shakeTolerance
= 0.0001;
346 gmx_bool shakeUseSOR
= false;
349 int lincslincsExpansionOrder
= 4;
350 real lincsWarnAngle
= 30.0;
352 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
353 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
354 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
355 lincslincsExpansionOrder
, lincsWarnAngle
);
357 std::string algorithmName
;
358 std::tie(pbcName
, algorithmName
) = GetParam();
359 t_pbc pbc
= pbcs_
.at(pbcName
);
362 algorithms_
.at(algorithmName
)(testData
.get(), pbc
);
364 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
365 checkConstrainsDirection(*testData
, pbc
);
366 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
367 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
369 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
372 TEST_P(ConstraintsTest
, TwoDisjointConstraints
)
375 std::string title
= "two disjoint constraints";
377 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
378 std::vector
<int> constraints
= { 0, 0, 1, 1, 2, 3 };
379 std::vector
<real
> constraintsR0
= { 2.0, 1.0 };
382 std::vector
<RVec
> x
= {
383 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
386 std::vector
<RVec
> xPrime
= {
387 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
390 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 } };
392 tensor virialScaledRef
= { { 3.3e-03, -1.7e-04, 5.6e-04 },
393 { -1.7e-04, 8.9e-06, -2.8e-05 },
394 { 5.6e-04, -2.8e-05, 8.9e-05 } };
396 real shakeTolerance
= 0.0001;
397 gmx_bool shakeUseSOR
= false;
400 int lincslincsExpansionOrder
= 4;
401 real lincsWarnAngle
= 30.0;
403 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
404 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
405 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
406 lincslincsExpansionOrder
, lincsWarnAngle
);
409 std::string algorithmName
;
410 std::tie(pbcName
, algorithmName
) = GetParam();
411 t_pbc pbc
= pbcs_
.at(pbcName
);
414 algorithms_
.at(algorithmName
)(testData
.get(), pbc
);
416 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
417 checkConstrainsDirection(*testData
, pbc
);
418 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
419 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
421 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
424 TEST_P(ConstraintsTest
, ThreeSequentialConstraints
)
427 std::string title
= "three atoms, connected longitudinally (e.g. CH2)";
429 std::vector
<real
> masses
= { 1.0, 12.0, 16.0 };
430 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2 };
431 std::vector
<real
> constraintsR0
= { 0.1, 0.2 };
433 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
434 real twoTenthsOverSqrtThree
= 0.2_real
/ std::sqrt(3.0_real
);
436 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, oneTenthOverSqrtTwo
, 0.0 },
438 { twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
, twoTenthsOverSqrtThree
} };
440 std::vector
<RVec
> xPrime
= { { 0.08, 0.07, 0.01 }, { -0.02, 0.01, -0.02 }, { 0.10, 0.12, 0.11 } };
442 std::vector
<RVec
> v
= { { 1.0, 0.0, 0.0 }, { 0.0, 1.0, 0.0 }, { 0.0, 0.0, 1.0 } };
444 tensor virialScaledRef
= { { 4.14e-03, 4.14e-03, 3.31e-03 },
445 { 4.14e-03, 4.14e-03, 3.31e-03 },
446 { 3.31e-03, 3.31e-03, 3.31e-03 } };
448 real shakeTolerance
= 0.0001;
449 gmx_bool shakeUseSOR
= false;
452 int lincslincsExpansionOrder
= 4;
453 real lincsWarnAngle
= 30.0;
455 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
456 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
457 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
458 lincslincsExpansionOrder
, lincsWarnAngle
);
461 std::string algorithmName
;
462 std::tie(pbcName
, algorithmName
) = GetParam();
463 t_pbc pbc
= pbcs_
.at(pbcName
);
466 algorithms_
.at(algorithmName
)(testData
.get(), pbc
);
468 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
469 checkConstrainsDirection(*testData
, pbc
);
470 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
471 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
473 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
476 TEST_P(ConstraintsTest
, ThreeConstraintsWithCentralAtom
)
479 std::string title
= "three atoms, connected to the central atom (e.g. CH3)";
481 std::vector
<real
> masses
= { 12.0, 1.0, 1.0, 1.0 };
482 std::vector
<int> constraints
= { 0, 0, 1, 0, 0, 2, 0, 0, 3 };
483 std::vector
<real
> constraintsR0
= { 0.1 };
486 std::vector
<RVec
> x
= {
487 { 0.00, 0.00, 0.00 }, { 0.10, 0.00, 0.00 }, { 0.00, -0.10, 0.00 }, { 0.00, 0.00, 0.10 }
490 std::vector
<RVec
> xPrime
= { { 0.004, 0.009, -0.010 },
491 { 0.110, -0.006, 0.003 },
492 { -0.007, -0.102, -0.007 },
493 { -0.005, 0.011, 0.102 } };
495 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 } };
497 tensor virialScaledRef
= { { 7.14e-04, 0.00e+00, 0.00e+00 },
498 { 0.00e+00, 1.08e-03, 0.00e+00 },
499 { 0.00e+00, 0.00e+00, 1.15e-03 } };
501 real shakeTolerance
= 0.0001;
502 gmx_bool shakeUseSOR
= false;
505 int lincslincsExpansionOrder
= 4;
506 real lincsWarnAngle
= 30.0;
508 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
509 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
510 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
511 lincslincsExpansionOrder
, lincsWarnAngle
);
514 std::string algorithmName
;
515 std::tie(pbcName
, algorithmName
) = GetParam();
516 t_pbc pbc
= pbcs_
.at(pbcName
);
519 algorithms_
.at(algorithmName
)(testData
.get(), pbc
);
521 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
522 checkConstrainsDirection(*testData
, pbc
);
523 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
524 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
526 checkVirialTensor(absoluteTolerance(0.0001), *testData
);
529 TEST_P(ConstraintsTest
, FourSequentialConstraints
)
532 std::string title
= "four atoms, connected longitudinally";
534 std::vector
<real
> masses
= { 0.5, 1.0 / 3.0, 0.25, 1.0 };
535 std::vector
<int> constraints
= { 0, 0, 1, 1, 1, 2, 2, 2, 3 };
536 std::vector
<real
> constraintsR0
= { 2.0, 1.0, 1.0 };
539 std::vector
<RVec
> x
= {
540 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
543 std::vector
<RVec
> xPrime
= {
544 { 2.50, -3.10, 15.70 }, { 0.51, -3.02, 15.55 }, { -0.50, -3.00, 15.20 }, { -1.51, -2.95, 15.05 }
547 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 } };
549 tensor virialScaledRef
= { { 1.15e-01, -4.20e-03, 2.12e-02 },
550 { -4.20e-03, 1.70e-04, -6.41e-04 },
551 { 2.12e-02, -6.41e-04, 5.45e-03 } };
553 real shakeTolerance
= 0.0001;
554 gmx_bool shakeUseSOR
= false;
557 int lincslincsExpansionOrder
= 8;
558 real lincsWarnAngle
= 30.0;
560 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
561 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
562 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
563 lincslincsExpansionOrder
, lincsWarnAngle
);
566 std::string algorithmName
;
567 std::tie(pbcName
, algorithmName
) = GetParam();
568 t_pbc pbc
= pbcs_
.at(pbcName
);
571 algorithms_
.at(algorithmName
)(testData
.get(), pbc
);
573 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
574 checkConstrainsDirection(*testData
, pbc
);
575 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
576 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
578 checkVirialTensor(absoluteTolerance(0.01), *testData
);
581 TEST_P(ConstraintsTest
, TriangleOfConstraints
)
584 std::string title
= "basic triangle (tree atoms, connected to each other)";
586 std::vector
<real
> masses
= { 1.0, 1.0, 1.0 };
587 std::vector
<int> constraints
= { 0, 0, 1, 2, 0, 2, 1, 1, 2 };
588 std::vector
<real
> constraintsR0
= { 0.1, 0.1, 0.1 };
590 real oneTenthOverSqrtTwo
= 0.1_real
/ std::sqrt(2.0_real
);
592 std::vector
<RVec
> x
= { { oneTenthOverSqrtTwo
, 0.0, 0.0 },
593 { 0.0, oneTenthOverSqrtTwo
, 0.0 },
594 { 0.0, 0.0, oneTenthOverSqrtTwo
} };
596 std::vector
<RVec
> xPrime
= { { 0.09, -0.02, 0.01 }, { -0.02, 0.10, -0.02 }, { 0.03, -0.01, 0.07 } };
598 std::vector
<RVec
> v
= { { 1.0, 1.0, 1.0 }, { -2.0, -2.0, -2.0 }, { 1.0, 1.0, 1.0 } };
600 tensor virialScaledRef
= { { 6.00e-04, -1.61e-03, 1.01e-03 },
601 { -1.61e-03, 2.53e-03, -9.25e-04 },
602 { 1.01e-03, -9.25e-04, -8.05e-05 } };
604 real shakeTolerance
= 0.0001;
605 gmx_bool shakeUseSOR
= false;
608 int lincslincsExpansionOrder
= 4;
609 real lincsWarnAngle
= 30.0;
611 std::unique_ptr
<ConstraintsTestData
> testData
= std::make_unique
<ConstraintsTestData
>(
612 title
, numAtoms
, masses
, constraints
, constraintsR0
, true, virialScaledRef
, false, 0,
613 real(0.0), real(0.001), x
, xPrime
, v
, shakeTolerance
, shakeUseSOR
, lincsNIter
,
614 lincslincsExpansionOrder
, lincsWarnAngle
);
617 std::string runnerName
;
618 std::tie(pbcName
, runnerName
) = GetParam();
619 t_pbc pbc
= pbcs_
.at(pbcName
);
622 algorithms_
.at(runnerName
)(testData
.get(), pbc
);
624 checkConstrainsLength(absoluteTolerance(0.0002), *testData
, pbc
);
625 checkConstrainsDirection(*testData
, pbc
);
626 checkCOMCoordinates(absoluteTolerance(0.0001), *testData
);
627 checkCOMVelocity(absoluteTolerance(0.0001), *testData
);
629 checkVirialTensor(absoluteTolerance(0.00001), *testData
);
633 INSTANTIATE_TEST_CASE_P(WithParameters
,
635 ::testing::Combine(::testing::Values("PBCNone", "PBCXYZ"),
636 ::testing::ValuesIn(getRunnersNames())));