Change the file structure in constraints tests
[gromacs.git] / src / gromacs / mdlib / tests / constr.cpp
blob1d6c85affd032774325ddee296f419fdc6adfc2b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
35 /*! \internal \file
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
47 #include "gmxpre.h"
49 #include "config.h"
51 #include <assert.h>
53 #include <unordered_map>
54 #include <vector>
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"
68 namespace gmx
70 namespace test
72 namespace
75 /*! \brief The two-dimensional parameter space for test.
77 * The test will run for all possible combinations of accessible
78 * values of the:
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");
96 return runnersNames;
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
115 * pre-computed data.
117 class ConstraintsTest : public ::testing::TestWithParam<ConstraintsTestParameters>
119 public:
120 //! PBC setups
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
137 t_pbc pbc;
139 // Infinitely small box
140 matrix boxNone = { {0, 0, 0}, {0, 0, 0}, {0, 0, 0} };
141 set_pbc(&pbc, epbcNONE, boxNone);
142 pbcs_["PBCNone"] = pbc;
144 // Rectangular box
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, epbcXYZ, boxXyz);
147 pbcs_["PBCXYZ"] = pbc;
150 // Algorithms
152 // SHAKE
153 algorithms_["SHAKE"] = applyShake;
154 // LINCS
155 algorithms_["LINCS"] = applyLincs;
156 // LINCS using CUDA (will only be called if CUDA is available)
157 algorithms_["LINCS_CUDA"] = applyLincsCuda;
160 /*! \brief
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 equal
164 * 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 (unsigned c = 0; c < testData.constraints_.size()/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);
179 RVec xij0, xij1;
180 real d0, d1;
181 if (pbc.ePBC == epbcXYZ)
183 pbc_dx_aiuc(&pbc, testData.x_[i], testData.x_[j], xij0);
184 pbc_dx_aiuc(&pbc, testData.xPrime_[i], testData.xPrime_[j], xij1);
186 else
188 rvec_sub(testData.x_[i], testData.x_[j], xij0);
189 rvec_sub(testData.xPrime_[i], testData.xPrime_[j], xij1);
191 d0 = norm(xij0);
192 d1 = norm(xij1);
193 EXPECT_REAL_EQ_TOL(r0, d1, tolerance) << gmx::formatString(
194 "rij = %f, which is not equal to r0 = %f for constraint #%u, between atoms %d and %d"
195 " (before constraining rij was %f).", d1, r0, c, i, j, d0);
199 /*! \brief
200 * The test on the final length of constrained bonds.
202 * Goes through all the constraints and checks if the direction of constraint has not changed
203 * by the algorithm (i.e. the constraints algorithm arrived to the solution that is closest
204 * to the initial system conformation).
206 * \param[in] testData Test data structure.
207 * \param[in] pbc Periodic boundary data.
209 void checkConstrainsDirection(const ConstraintsTestData &testData, t_pbc pbc)
212 for (unsigned c = 0; c < testData.constraints_.size()/3; c++)
214 int i = testData.constraints_.at(3*c + 1);
215 int j = testData.constraints_.at(3*c + 2);
216 RVec xij0, xij1;
217 if (pbc.ePBC == epbcXYZ)
219 pbc_dx_aiuc(&pbc, testData.x_[i], testData.x_[j], xij0);
220 pbc_dx_aiuc(&pbc, testData.xPrime_[i], testData.xPrime_[j], xij1);
222 else
224 rvec_sub(testData.x_[i], testData.x_[j], xij0);
225 rvec_sub(testData.xPrime_[i], testData.xPrime_[j], xij1);
228 real dot = xij0.dot(xij1);
230 EXPECT_GE(dot, 0.0) << gmx::formatString(
231 "The constraint %u changed direction. Constraining algorithm might have returned the wrong root "
232 "of the constraints equation.", c);
237 /*! \brief
238 * The test on the coordinates of the center of the mass (COM) of the system.
240 * Checks if the center of mass has not been shifted by the constraints. Note,
241 * that this test does not take into account the periodic boundary conditions.
242 * Hence it will not work should the constraints decide to move atoms across
243 * PBC borders.
245 * \param[in] tolerance Allowed tolerance in COM coordinates.
246 * \param[in] testData Test data structure.
248 void checkCOMCoordinates(FloatingPointTolerance tolerance, const ConstraintsTestData &testData)
251 RVec comPrime0({0.0, 0.0, 0.0});
252 RVec comPrime({0.0, 0.0, 0.0});
253 for (int i = 0; i < testData.numAtoms_; i++)
255 comPrime0 += testData.masses_[i]*testData.xPrime0_[i];
256 comPrime += testData.masses_[i]*testData.xPrime_[i];
259 comPrime0 /= testData.numAtoms_;
260 comPrime /= testData.numAtoms_;
262 EXPECT_REAL_EQ_TOL(comPrime[XX], comPrime0[XX], tolerance)
263 << "Center of mass was shifted by constraints in x-direction.";
264 EXPECT_REAL_EQ_TOL(comPrime[YY], comPrime0[YY], tolerance)
265 << "Center of mass was shifted by constraints in y-direction.";
266 EXPECT_REAL_EQ_TOL(comPrime[ZZ], comPrime0[ZZ], tolerance)
267 << "Center of mass was shifted by constraints in z-direction.";
271 /*! \brief
272 * The test on the velocity of the center of the mass (COM) of the system.
274 * Checks if the velocity of the center of mass has not changed.
276 * \param[in] tolerance Allowed tolerance in COM velocity components.
277 * \param[in] testData Test data structure.
279 void checkCOMVelocity(FloatingPointTolerance tolerance, const ConstraintsTestData &testData)
282 RVec comV0({0.0, 0.0, 0.0});
283 RVec comV({0.0, 0.0, 0.0});
284 for (int i = 0; i < testData.numAtoms_; i++)
286 comV0 += testData.masses_[i]*testData.v0_[i];
287 comV += testData.masses_[i]*testData.v_[i];
289 comV0 /= testData.numAtoms_;
290 comV /= testData.numAtoms_;
292 EXPECT_REAL_EQ_TOL(comV[XX], comV0[XX], tolerance)
293 << "Velocity of the center of mass in x-direction has been changed by constraints.";
294 EXPECT_REAL_EQ_TOL(comV[YY], comV0[YY], tolerance)
295 << "Velocity of the center of mass in y-direction has been changed by constraints.";
296 EXPECT_REAL_EQ_TOL(comV[ZZ], comV0[ZZ], tolerance)
297 << "Velocity of the center of mass in z-direction has been changed by constraints.";
300 /*! \brief
301 * The test of virial tensor.
303 * Checks if the values in the scaled virial tensor are equal to pre-computed values.
305 * \param[in] tolerance Tolerance for the tensor values.
306 * \param[in] testData Test data structure.
308 void checkVirialTensor(FloatingPointTolerance tolerance, const ConstraintsTestData &testData)
310 for (int i = 0; i < DIM; i++)
312 for (int j = 0; j < DIM; j++)
314 EXPECT_REAL_EQ_TOL(testData.virialScaledRef_[i][j], testData.virialScaled_[i][j],
315 tolerance) << gmx::formatString(
316 "Values in virial tensor at [%d][%d] are not within the tolerance from reference value.", i, j);
322 TEST_P(ConstraintsTest, SingleConstraint){
323 std::string title = "one constraint (e.g. OH)";
324 int numAtoms = 2;
326 std::vector<real> masses = {1.0, 12.0};
327 std::vector<int> constraints = {0, 0, 1};
328 std::vector<real> constraintsR0 = {0.1};
330 real oneTenthOverSqrtTwo = 0.1_real / std::sqrt(2.0_real);
332 std::vector<RVec> x = {{ 0.0, oneTenthOverSqrtTwo, 0.0 },
333 { oneTenthOverSqrtTwo, 0.0, 0.0 }};
334 std::vector<RVec> xPrime = {{ 0.01, 0.08, 0.01 },
335 { 0.06, 0.01, -0.01 }};
336 std::vector<RVec> v = {{ 1.0, 2.0, 3.0 },
337 { 3.0, 2.0, 1.0 }};
339 tensor virialScaledRef = {{-5.58e-04, 5.58e-04, 0.00e+00 },
340 { 5.58e-04, -5.58e-04, 0.00e+00 },
341 { 0.00e+00, 0.00e+00, 0.00e+00 }};
343 real shakeTolerance = 0.0001;
344 gmx_bool shakeUseSOR = false;
346 int lincsNIter = 1;
347 int lincslincsExpansionOrder = 4;
348 real lincsWarnAngle = 30.0;
350 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
351 (title, numAtoms, masses,
352 constraints, constraintsR0,
353 true, virialScaledRef,
354 false, 0,
355 real(0.0), real(0.001),
356 x, xPrime, v,
357 shakeTolerance, shakeUseSOR,
358 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
359 std::string pbcName;
360 std::string algorithmName;
361 std::tie(pbcName, algorithmName) = GetParam();
362 t_pbc pbc = pbcs_.at(pbcName);
364 // Apply constraints
365 algorithms_.at(algorithmName)(testData.get(), pbc);
367 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
368 checkConstrainsDirection(*testData, pbc);
369 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
370 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
372 checkVirialTensor(absoluteTolerance(0.0001), *testData);
376 TEST_P(ConstraintsTest, TwoDisjointConstraints){
378 std::string title = "two disjoint constraints";
379 int numAtoms = 4;
380 std::vector<real> masses = {0.5, 1.0/3.0, 0.25, 1.0};
381 std::vector<int> constraints = {0, 0, 1, 1, 2, 3};
382 std::vector<real> constraintsR0 = {2.0, 1.0};
385 std::vector<RVec> x = {{ 2.50, -3.10, 15.70 },
386 { 0.51, -3.02, 15.55 },
387 { -0.50, -3.00, 15.20 },
388 { -1.51, -2.95, 15.05 }};
390 std::vector<RVec> xPrime = {{ 2.50, -3.10, 15.70 },
391 { 0.51, -3.02, 15.55 },
392 { -0.50, -3.00, 15.20 },
393 { -1.51, -2.95, 15.05 }};
395 std::vector<RVec> v = {{ 0.0, 1.0, 0.0 },
396 { 1.0, 0.0, 0.0 },
397 { 0.0, 0.0, 1.0 },
398 { 0.0, 0.0, 0.0 }};
400 tensor virialScaledRef = {{ 3.3e-03, -1.7e-04, 5.6e-04 },
401 {-1.7e-04, 8.9e-06, -2.8e-05 },
402 { 5.6e-04, -2.8e-05, 8.9e-05 }};
404 real shakeTolerance = 0.0001;
405 gmx_bool shakeUseSOR = false;
407 int lincsNIter = 1;
408 int lincslincsExpansionOrder = 4;
409 real lincsWarnAngle = 30.0;
411 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
412 (title, numAtoms, masses,
413 constraints, constraintsR0,
414 true, virialScaledRef,
415 false, 0,
416 real(0.0), real(0.001),
417 x, xPrime, v,
418 shakeTolerance, shakeUseSOR,
419 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
421 std::string pbcName;
422 std::string algorithmName;
423 std::tie(pbcName, algorithmName) = GetParam();
424 t_pbc pbc = pbcs_.at(pbcName);
426 // Apply constraints
427 algorithms_.at(algorithmName)(testData.get(), pbc);
429 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
430 checkConstrainsDirection(*testData, pbc);
431 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
432 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
434 checkVirialTensor(absoluteTolerance(0.0001), *testData);
438 TEST_P(ConstraintsTest, ThreeSequentialConstraints){
440 std::string title = "three atoms, connected longitudinally (e.g. CH2)";
441 int numAtoms = 3;
442 std::vector<real> masses = {1.0, 12.0, 16.0 };
443 std::vector<int> constraints = {0, 0, 1, 1, 1, 2};
444 std::vector<real> constraintsR0 = {0.1, 0.2};
446 real oneTenthOverSqrtTwo = 0.1_real / std::sqrt(2.0_real);
447 real twoTenthsOverSqrtThree = 0.2_real / std::sqrt(3.0_real);
449 std::vector<RVec> x = {{ oneTenthOverSqrtTwo, oneTenthOverSqrtTwo, 0.0 },
450 { 0.0, 0.0, 0.0 },
451 { twoTenthsOverSqrtThree, twoTenthsOverSqrtThree, twoTenthsOverSqrtThree }};
453 std::vector<RVec> xPrime = {{ 0.08, 0.07, 0.01 },
454 { -0.02, 0.01, -0.02 },
455 { 0.10, 0.12, 0.11 }};
457 std::vector<RVec> v = {{ 1.0, 0.0, 0.0 },
458 { 0.0, 1.0, 0.0 },
459 { 0.0, 0.0, 1.0 }};
461 tensor virialScaledRef = {{ 4.14e-03, 4.14e-03, 3.31e-03},
462 { 4.14e-03, 4.14e-03, 3.31e-03},
463 { 3.31e-03, 3.31e-03, 3.31e-03}};
465 real shakeTolerance = 0.0001;
466 gmx_bool shakeUseSOR = false;
468 int lincsNIter = 1;
469 int lincslincsExpansionOrder = 4;
470 real lincsWarnAngle = 30.0;
472 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
473 (title, numAtoms, masses,
474 constraints, constraintsR0,
475 true, virialScaledRef,
476 false, 0,
477 real(0.0), real(0.001),
478 x, xPrime, v,
479 shakeTolerance, shakeUseSOR,
480 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
482 std::string pbcName;
483 std::string algorithmName;
484 std::tie(pbcName, algorithmName) = GetParam();
485 t_pbc pbc = pbcs_.at(pbcName);
487 // Apply constraints
488 algorithms_.at(algorithmName)(testData.get(), pbc);
490 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
491 checkConstrainsDirection(*testData, pbc);
492 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
493 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
495 checkVirialTensor(absoluteTolerance(0.0001), *testData);
499 TEST_P(ConstraintsTest, ThreeConstraintsWithCentralAtom){
501 std::string title = "three atoms, connected to the central atom (e.g. CH3)";
502 int numAtoms = 4;
503 std::vector<real> masses = {12.0, 1.0, 1.0, 1.0 };
504 std::vector<int> constraints = {0, 0, 1, 0, 0, 2, 0, 0, 3};
505 std::vector<real> constraintsR0 = {0.1};
508 std::vector<RVec> x = {{ 0.00, 0.00, 0.00 },
509 { 0.10, 0.00, 0.00 },
510 { 0.00, -0.10, 0.00 },
511 { 0.00, 0.00, 0.10 }};
513 std::vector<RVec> xPrime = {{ 0.004, 0.009, -0.010 },
514 { 0.110, -0.006, 0.003 },
515 {-0.007, -0.102, -0.007 },
516 {-0.005, 0.011, 0.102 }};
518 std::vector<RVec> v = {{ 1.0, 0.0, 0.0 },
519 { 1.0, 0.0, 0.0 },
520 { 1.0, 0.0, 0.0 },
521 { 1.0, 0.0, 0.0 }};
523 tensor virialScaledRef = {{7.14e-04, 0.00e+00, 0.00e+00},
524 {0.00e+00, 1.08e-03, 0.00e+00},
525 {0.00e+00, 0.00e+00, 1.15e-03}};
527 real shakeTolerance = 0.0001;
528 gmx_bool shakeUseSOR = false;
530 int lincsNIter = 1;
531 int lincslincsExpansionOrder = 4;
532 real lincsWarnAngle = 30.0;
534 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
535 (title, numAtoms, masses,
536 constraints, constraintsR0,
537 true, virialScaledRef,
538 false, 0,
539 real(0.0), real(0.001),
540 x, xPrime, v,
541 shakeTolerance, shakeUseSOR,
542 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
544 std::string pbcName;
545 std::string algorithmName;
546 std::tie(pbcName, algorithmName) = GetParam();
547 t_pbc pbc = pbcs_.at(pbcName);
549 // Apply constraints
550 algorithms_.at(algorithmName)(testData.get(), pbc);
552 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
553 checkConstrainsDirection(*testData, pbc);
554 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
555 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
557 checkVirialTensor(absoluteTolerance(0.0001), *testData);
560 TEST_P(ConstraintsTest, FourSequentialConstraints){
562 std::string title = "four atoms, connected longitudinally";
563 int numAtoms = 4;
564 std::vector<real> masses = {0.5, 1.0/3.0, 0.25, 1.0};
565 std::vector<int> constraints = {0, 0, 1, 1, 1, 2, 2, 2, 3};
566 std::vector<real> constraintsR0 = {2.0, 1.0, 1.0};
569 std::vector<RVec> x = {{ 2.50, -3.10, 15.70 },
570 { 0.51, -3.02, 15.55 },
571 { -0.50, -3.00, 15.20 },
572 { -1.51, -2.95, 15.05 }};
574 std::vector<RVec> xPrime = {{ 2.50, -3.10, 15.70 },
575 { 0.51, -3.02, 15.55 },
576 { -0.50, -3.00, 15.20 },
577 { -1.51, -2.95, 15.05 }};
579 std::vector<RVec> v = {{ 0.0, 0.0, 2.0 },
580 { 0.0, 0.0, 3.0 },
581 { 0.0, 0.0, -4.0 },
582 { 0.0, 0.0, -1.0 }};
584 tensor virialScaledRef = {{ 1.15e-01, -4.20e-03, 2.12e-02},
585 {-4.20e-03, 1.70e-04, -6.41e-04},
586 { 2.12e-02, -6.41e-04, 5.45e-03}};
588 real shakeTolerance = 0.0001;
589 gmx_bool shakeUseSOR = false;
591 int lincsNIter = 4;
592 int lincslincsExpansionOrder = 8;
593 real lincsWarnAngle = 30.0;
595 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
596 (title, numAtoms, masses,
597 constraints, constraintsR0,
598 true, virialScaledRef,
599 false, 0,
600 real(0.0), real(0.001),
601 x, xPrime, v,
602 shakeTolerance, shakeUseSOR,
603 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
605 std::string pbcName;
606 std::string algorithmName;
607 std::tie(pbcName, algorithmName) = GetParam();
608 t_pbc pbc = pbcs_.at(pbcName);
610 // Apply constraints
611 algorithms_.at(algorithmName)(testData.get(), pbc);
613 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
614 checkConstrainsDirection(*testData, pbc);
615 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
616 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
618 checkVirialTensor(absoluteTolerance(0.01), *testData);
622 TEST_P(ConstraintsTest, TriangleOfConstraints){
624 std::string title = "basic triangle (tree atoms, connected to each other)";
625 int numAtoms = 3;
626 std::vector<real> masses = {1.0, 1.0, 1.0};
627 std::vector<int> constraints = {0, 0, 1, 2, 0, 2, 1, 1, 2};
628 std::vector<real> constraintsR0 = {0.1, 0.1, 0.1};
630 real oneTenthOverSqrtTwo = 0.1_real / std::sqrt(2.0_real);
632 std::vector<RVec> x = {{ oneTenthOverSqrtTwo, 0.0, 0.0 },
633 { 0.0, oneTenthOverSqrtTwo, 0.0 },
634 { 0.0, 0.0, oneTenthOverSqrtTwo }};
636 std::vector<RVec> xPrime = {{ 0.09, -0.02, 0.01 },
637 { -0.02, 0.10, -0.02 },
638 { 0.03, -0.01, 0.07 }};
640 std::vector<RVec> v = {{ 1.0, 1.0, 1.0 },
641 { -2.0, -2.0, -2.0 },
642 { 1.0, 1.0, 1.0 }};
644 tensor virialScaledRef = {{ 6.00e-04, -1.61e-03, 1.01e-03},
645 {-1.61e-03, 2.53e-03, -9.25e-04},
646 { 1.01e-03, -9.25e-04, -8.05e-05}};
648 real shakeTolerance = 0.0001;
649 gmx_bool shakeUseSOR = false;
651 int lincsNIter = 1;
652 int lincslincsExpansionOrder = 4;
653 real lincsWarnAngle = 30.0;
655 std::unique_ptr<ConstraintsTestData> testData = std::make_unique<ConstraintsTestData>
656 (title, numAtoms, masses,
657 constraints, constraintsR0,
658 true, virialScaledRef,
659 false, 0,
660 real(0.0), real(0.001),
661 x, xPrime, v,
662 shakeTolerance, shakeUseSOR,
663 lincsNIter, lincslincsExpansionOrder, lincsWarnAngle);
665 std::string pbcName;
666 std::string runnerName;
667 std::tie(pbcName, runnerName) = GetParam();
668 t_pbc pbc = pbcs_.at(pbcName);
670 // Apply constraints
671 algorithms_.at(runnerName)(testData.get(), pbc);
673 checkConstrainsLength(absoluteTolerance(0.0002), *testData, pbc);
674 checkConstrainsDirection(*testData, pbc);
675 checkCOMCoordinates(absoluteTolerance(0.0001), *testData);
676 checkCOMVelocity(absoluteTolerance(0.0001), *testData);
678 checkVirialTensor(absoluteTolerance(0.00001), *testData);
683 INSTANTIATE_TEST_CASE_P(WithParameters, ConstraintsTest,
684 ::testing::Combine(::testing::Values("PBCNone", "PBCXYZ"),
685 ::testing::ValuesIn(getRunnersNames())));
687 } // namespace
688 } // namespace test
689 } // namespace gmx