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.
37 * Implements test of bonded force routines
39 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
40 * \ingroup module_listed_forces
44 #include "gromacs/listed_forces/bonded.h"
49 #include <unordered_map>
51 #include <gtest/gtest.h>
53 #include "gromacs/listed_forces/listed_forces.h"
54 #include "gromacs/math/units.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/mdtypes/mdatom.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/topology/idef.h"
61 #include "gromacs/utility/strconvert.h"
63 #include "testutils/refdata.h"
64 #include "testutils/testasserts.h"
71 //! Number of atoms used in these tests.
72 constexpr int c_numAtoms = 4;
74 /*! \brief Output from bonded kernels
76 * \todo Later this might turn into the actual output struct. */
77 struct OutputQuantities
79 //! Energy of this interaction
81 //! Derivative with respect to lambda
84 rvec fshift[N_IVEC] = { { 0 } };
86 rvec4 f[c_numAtoms] = { { 0 } };
89 /*! \brief Utility to check the output from bonded tests
91 * \param[in] checker Reference checker
92 * \param[in] output The output from the test to check
94 void checkOutput(test::TestReferenceChecker* checker, const OutputQuantities& output)
96 checker->checkReal(output.energy, "Epot ");
97 // Should still be zero if not doing FEP, so may as well test it.
98 checker->checkReal(output.dvdlambda, "dVdlambda ");
99 checker->checkVector(output.fshift[CENTRAL], "Central shift forces");
100 checker->checkSequence(std::begin(output.f), std::end(output.f), "Forces");
103 /*! \brief Input structure for listed forces tests
110 //! Tolerance for float evaluation
112 //! Tolerance for double evaluation
113 double dtoler = 1e-8;
114 //! Do free energy perturbation?
116 //! Interaction parameters
117 t_iparams iparams = { { 0 } };
122 /*! \brief Constructor with tolerance
124 * \param[in] ftol Single precision tolerance
125 * \param[in] dtol Double precision tolerance
127 iListInput(float ftol, double dtol)
132 /*! \brief Set parameters for harmonic potential
134 * Free energy perturbation is turned on when A
135 * and B parameters are different.
136 * \param[in] ft Function type
137 * \param[in] rA Equilibrium value A
138 * \param[in] krA Force constant A
139 * \param[in] rB Equilibrium value B
140 * \param[in] krB Force constant B
141 * \return The structure itself.
143 iListInput setHarmonic(int ft, real rA, real krA, real rB, real krB)
145 iparams.harmonic.rA = rA;
146 iparams.harmonic.rB = rB;
147 iparams.harmonic.krA = krA;
148 iparams.harmonic.krB = krB;
150 fep = (rA != rB || krA != krB);
153 /*! \brief Set parameters for harmonic potential
155 * \param[in] ft Function type
156 * \param[in] rA Equilibrium value
157 * \param[in] krA Force constant
158 * \return The structure itself.
160 iListInput setHarmonic(int ft, real rA, real krA) { return setHarmonic(ft, rA, krA, rA, krA); }
161 /*! \brief Set parameters for cubic potential
163 * \param[in] b0 Equilibrium bond length
164 * \param[in] kb Harmonic force constant
165 * \param[in] kcub Cubic force constant
166 * \return The structure itself.
168 iListInput setCubic(real b0, real kb, real kcub)
170 ftype = F_CUBICBONDS;
171 iparams.cubic.b0 = b0;
172 iparams.cubic.kb = kb;
173 iparams.cubic.kcub = kcub;
176 /*! \brief Set parameters for morse potential
178 * Free energy perturbation is turned on when A
179 * and B parameters are different.
180 * \param[in] b0A Equilibrium value A
181 * \param[in] cbA Force constant A
182 * \param[in] betaA Steepness parameter A
183 * \param[in] b0B Equilibrium value B
184 * \param[in] cbB Force constant B
185 * \param[in] betaB Steepness parameter B
186 * \return The structure itself.
188 iListInput setMorse(real b0A, real cbA, real betaA, real b0B, real cbB, real betaB)
191 iparams.morse.b0A = b0A;
192 iparams.morse.cbA = cbA;
193 iparams.morse.betaA = betaA;
194 iparams.morse.b0B = b0B;
195 iparams.morse.cbB = cbB;
196 iparams.morse.betaB = betaB;
197 fep = (b0A != b0B || cbA != cbB || betaA != betaB);
200 /*! \brief Set parameters for morse potential
202 * \param[in] b0A Equilibrium value
203 * \param[in] cbA Force constant
204 * \param[in] betaA Steepness parameter
205 * \return The structure itself.
207 iListInput setMorse(real b0A, real cbA, real betaA)
209 return setMorse(b0A, cbA, betaA, b0A, cbA, betaA);
211 /*! \brief Set parameters for fene potential
213 * \param[in] bm Equilibrium bond length
214 * \param[in] kb Force constant
215 * \return The structure itself.
217 iListInput setFene(real bm, real kb)
220 iparams.fene.bm = bm;
221 iparams.fene.kb = kb;
224 /*! \brief Set parameters for linear angle potential
226 * Free energy perturbation is turned on when A
227 * and B parameters are different.
228 * \param[in] klinA Force constant A
229 * \param[in] aA The position of the central atom A
230 * \param[in] klinB Force constant B
231 * \param[in] aB The position of the central atom B
232 * \return The structure itself.
234 iListInput setLinearAngle(real klinA, real aA, real klinB, real aB)
236 ftype = F_LINEAR_ANGLES;
237 iparams.linangle.klinA = klinA;
238 iparams.linangle.aA = aA;
239 iparams.linangle.klinB = klinB;
240 iparams.linangle.aB = aB;
241 fep = (klinA != klinB || aA != aB);
244 /*! \brief Set parameters for linear angle potential
246 * \param[in] klinA Force constant
247 * \param[in] aA The position of the central atom
248 * \return The structure itself.
250 iListInput setLinearAngle(real klinA, real aA) { return setLinearAngle(klinA, aA, klinA, aA); }
251 /*! \brief Set parameters for Urey Bradley potential
253 * Free energy perturbation is turned on when A
254 * and B parameters are different.
255 * \param[in] thetaA Equilibrium angle A
256 * \param[in] kthetaA Force constant A
257 * \param[in] r13A The distance between i and k atoms A
258 * \param[in] kUBA The force constant for 1-3 distance A
259 * \param[in] thetaB Equilibrium angle B
260 * \param[in] kthetaB Force constant B
261 * \param[in] r13B The distance between i and k atoms B
262 * \param[in] kUBB The force constant for 1-3 distance B
263 * \return The structure itself.
266 setUreyBradley(real thetaA, real kthetaA, real r13A, real kUBA, real thetaB, real kthetaB, real r13B, real kUBB)
268 ftype = F_UREY_BRADLEY;
269 iparams.u_b.thetaA = thetaA;
270 iparams.u_b.kthetaA = kthetaA;
271 iparams.u_b.r13A = r13A;
272 iparams.u_b.kUBA = kUBA;
273 iparams.u_b.thetaB = thetaB;
274 iparams.u_b.kthetaB = kthetaB;
275 iparams.u_b.r13B = r13B;
276 iparams.u_b.kUBB = kUBB;
277 fep = (thetaA != thetaB || kthetaA != kthetaB || r13A != r13B || kUBA != kUBB);
280 /*! \brief Set parameters for Urey Bradley potential
282 * \param[in] thetaA Equilibrium angle
283 * \param[in] kthetaA Force constant
284 * \param[in] r13A The distance between i and k atoms
285 * \param[in] kUBA The force constant for 1-3 distance
286 * \return The structure itself.
288 iListInput setUreyBradley(real thetaA, real kthetaA, real r13A, real kUBA)
290 return setUreyBradley(thetaA, kthetaA, r13A, kUBA, thetaA, kthetaA, r13A, kUBA);
292 /*! \brief Set parameters for Cross Bond Bonds potential
294 * \param[in] r1e First bond length i-j
295 * \param[in] r2e Second bond length i-k
296 * \param[in] krr The force constant
297 * \return The structure itself.
299 iListInput setCrossBondBonds(real r1e, real r2e, real krr)
301 ftype = F_CROSS_BOND_BONDS;
302 iparams.cross_bb.r1e = r1e;
303 iparams.cross_bb.r2e = r2e;
304 iparams.cross_bb.krr = krr;
307 /*! \brief Set parameters for Cross Bond Angles potential
309 * \param[in] r1e First bond length i-j
310 * \param[in] r2e Second bond length j-k
311 * \param[in] r3e Third bond length i-k
312 * \param[in] krt The force constant
313 * \return The structure itself.
315 iListInput setCrossBondAngles(real r1e, real r2e, real r3e, real krt)
317 ftype = F_CROSS_BOND_ANGLES;
318 iparams.cross_ba.r1e = r1e;
319 iparams.cross_ba.r2e = r2e;
320 iparams.cross_ba.r3e = r3e;
321 iparams.cross_ba.krt = krt;
324 /*! \brief Set parameters for Quartic Angles potential
326 * \param[in] theta Angle
327 * \param[in] c Array of parameters
328 * \return The structure itself.
330 iListInput setQuarticAngles(real theta, const real c[5])
332 ftype = F_QUARTIC_ANGLES;
333 iparams.qangle.theta = theta;
334 iparams.qangle.c[0] = c[0];
335 iparams.qangle.c[1] = c[1];
336 iparams.qangle.c[2] = c[2];
337 iparams.qangle.c[3] = c[3];
338 iparams.qangle.c[4] = c[4];
341 /*! \brief Set parameters for proper dihedrals potential
343 * Free energy perturbation is turned on when A
344 * and B parameters are different.
345 * \param[in] ft Function type
346 * \param[in] phiA Dihedral angle A
347 * \param[in] cpA Force constant A
348 * \param[in] mult Multiplicity of the angle
349 * \param[in] phiB Dihedral angle B
350 * \param[in] cpB Force constant B
351 * \return The structure itself.
353 iListInput setPDihedrals(int ft, real phiA, real cpA, int mult, real phiB, real cpB)
356 iparams.pdihs.phiA = phiA;
357 iparams.pdihs.cpA = cpA;
358 iparams.pdihs.phiB = phiB;
359 iparams.pdihs.cpB = cpB;
360 iparams.pdihs.mult = mult;
361 fep = (phiA != phiB || cpA != cpB);
364 /*! \brief Set parameters for proper dihedrals potential
366 * \param[in] ft Function type
367 * \param[in] phiA Dihedral angle
368 * \param[in] cpA Force constant
369 * \param[in] mult Multiplicity of the angle
370 * \return The structure itself.
372 iListInput setPDihedrals(int ft, real phiA, real cpA, int mult)
374 return setPDihedrals(ft, phiA, cpA, mult, phiA, cpA);
376 /*! \brief Set parameters for Ryckaert-Bellemans dihedrals potential
378 * Free energy perturbation is turned on when A
379 * and B parameters are different.
380 * \param[in] rbcA Force constants A
381 * \param[in] rbcB Force constants B
382 * \return The structure itself.
384 iListInput setRbDihedrals(const real rbcA[NR_RBDIHS], const real rbcB[NR_RBDIHS])
388 for (int i = 0; i < NR_RBDIHS; i++)
390 iparams.rbdihs.rbcA[i] = rbcA[i];
391 iparams.rbdihs.rbcB[i] = rbcB[i];
392 fep = fep || (rbcA[i] != rbcB[i]);
396 /*! \brief Set parameters for Ryckaert-Bellemans dihedrals potential
398 * \param[in] rbc Force constants
399 * \return The structure itself.
401 iListInput setRbDihedrals(const real rbc[NR_RBDIHS]) { return setRbDihedrals(rbc, rbc); }
402 /*! \brief Set parameters for Polarization
404 * \param[in] alpha Polarizability
405 * \return The structure itself.
407 iListInput setPolarization(real alpha)
409 ftype = F_POLARIZATION;
411 iparams.polarize.alpha = alpha;
414 /*! \brief Set parameters for Anharmonic Polarization
416 * \param[in] alpha Polarizability (nm^3)
417 * \param[in] drcut The cut-off distance (nm) after which the
418 * fourth power kicks in
419 * \param[in] khyp The force constant for the fourth power
420 * \return The structure itself.
422 iListInput setAnharmPolarization(real alpha, real drcut, real khyp)
424 ftype = F_ANHARM_POL;
426 iparams.anharm_polarize.alpha = alpha;
427 iparams.anharm_polarize.drcut = drcut;
428 iparams.anharm_polarize.khyp = khyp;
431 /*! \brief Set parameters for Thole Polarization
433 * \param[in] a Thole factor
434 * \param[in] alpha1 Polarizability 1 (nm^3)
435 * \param[in] alpha2 Polarizability 2 (nm^3)
436 * \param[in] rfac Distance factor
437 * \return The structure itself.
439 iListInput setTholePolarization(real a, real alpha1, real alpha2, real rfac)
444 iparams.thole.alpha1 = alpha1;
445 iparams.thole.alpha2 = alpha2;
446 iparams.thole.rfac = rfac;
449 /*! \brief Set parameters for Water Polarization
451 * \param[in] alpha_x Polarizability X (nm^3)
452 * \param[in] alpha_y Polarizability Y (nm^3)
453 * \param[in] alpha_z Polarizability Z (nm^3)
454 * \param[in] rOH Oxygen-Hydrogen distance
455 * \param[in] rHH Hydrogen-Hydrogen distance
456 * \param[in] rOD Oxygen-Dummy distance
457 * \return The structure itself.
459 iListInput setWaterPolarization(real alpha_x, real alpha_y, real alpha_z, real rOH, real rHH, real rOD)
463 iparams.wpol.al_x = alpha_x;
464 iparams.wpol.al_y = alpha_y;
465 iparams.wpol.al_z = alpha_z;
466 iparams.wpol.rOH = rOH;
467 iparams.wpol.rHH = rHH;
468 iparams.wpol.rOD = rOD;
473 /*! \brief Utility to fill iatoms struct
475 * \param[in] ftype Function type
476 * \param[out] iatoms Pointer to iatoms struct
478 void fillIatoms(int ftype, std::vector<t_iatom>* iatoms)
480 std::unordered_map<int, std::vector<int>> ia = { { 2, { 0, 0, 1, 0, 1, 2, 0, 2, 3 } },
481 { 3, { 0, 0, 1, 2, 0, 1, 2, 3 } },
482 { 4, { 0, 0, 1, 2, 3 } },
483 { 5, { 0, 0, 1, 2, 3, 0 } } };
484 EXPECT_TRUE(ftype >= 0 && ftype < F_NRE);
485 int nral = interaction_function[ftype].nratoms;
486 for (auto& i : ia[nral])
488 iatoms->push_back(i);
492 class ListedForcesTest :
493 public ::testing::TestWithParam<std::tuple<iListInput, std::vector<gmx::RVec>, int>>
498 std::vector<gmx::RVec> x_;
501 test::TestReferenceData refData_;
502 test::TestReferenceChecker checker_;
503 ListedForcesTest() : checker_(refData_.rootChecker())
505 input_ = std::get<0>(GetParam());
506 x_ = std::get<1>(GetParam());
507 epbc_ = std::get<2>(GetParam());
509 box_[0][0] = box_[1][1] = box_[2][2] = 1.5;
510 set_pbc(&pbc_, epbc_, box_);
511 // We need quite specific tolerances here since angle functions
512 // etc. are not very precise and reproducible.
513 test::FloatingPointTolerance tolerance(test::FloatingPointTolerance(
514 input_.ftoler, 1.0e-6, input_.dtoler, 1.0e-12, 10000, 100, false));
515 checker_.setDefaultTolerance(tolerance);
517 void testOneIfunc(test::TestReferenceChecker* checker, const std::vector<t_iatom>& iatoms, const real lambda)
519 SCOPED_TRACE(std::string("Testing PBC ") + epbc_names[epbc_]);
520 std::vector<int> ddgatindex = { 0, 1, 2, 3 };
521 std::vector<real> chargeA = { 1.5, -2.0, 1.5, -1.0 };
522 t_mdatoms mdatoms = { 0 };
523 mdatoms.chargeA = chargeA.data();
524 OutputQuantities output;
525 output.energy = calculateSimpleBond(
526 input_.ftype, iatoms.size(), iatoms.data(), &input_.iparams,
527 as_rvec_array(x_.data()), output.f, output.fshift, &pbc_,
528 /* const struct t_graph *g */ nullptr, lambda, &output.dvdlambda, &mdatoms,
529 /* struct t_fcdata * */ nullptr, ddgatindex.data(),
530 BondedKernelFlavor::ForcesAndVirialAndEnergy);
531 // Internal consistency test of both test input
532 // and bonded functions.
533 EXPECT_TRUE((input_.fep || (output.dvdlambda == 0.0)));
534 checkOutput(checker, output);
538 test::TestReferenceChecker thisChecker =
539 checker_.checkCompound("FunctionType", interaction_function[input_.ftype].name)
540 .checkCompound("FEP", (input_.fep ? "Yes" : "No"));
541 std::vector<t_iatom> iatoms;
542 fillIatoms(input_.ftype, &iatoms);
545 const int numLambdas = 3;
546 for (int i = 0; i < numLambdas; ++i)
548 const real lambda = i / (numLambdas - 1.0);
549 auto valueChecker = thisChecker.checkCompound("Lambda", toString(lambda));
550 testOneIfunc(&valueChecker, iatoms, lambda);
555 testOneIfunc(&thisChecker, iatoms, 0.0);
560 TEST_P(ListedForcesTest, Ifunc)
565 //! Function types for testing bonds. Add new terms at the end.
566 std::vector<iListInput> c_InputBonds = {
567 { iListInput().setHarmonic(F_BONDS, 0.15, 500.0) },
568 { iListInput(2e-6F, 1e-8).setHarmonic(F_BONDS, 0.15, 500.0, 0.17, 400.0) },
569 { iListInput(1e-4F, 1e-8).setHarmonic(F_G96BONDS, 0.15, 50.0) },
570 { iListInput().setHarmonic(F_G96BONDS, 0.15, 50.0, 0.17, 40.0) },
571 { iListInput().setCubic(0.16, 50.0, 2.0) },
572 { iListInput(2e-6F, 1e-8).setMorse(0.15, 50.0, 2.0, 0.17, 40.0, 1.6) },
573 { iListInput(2e-6F, 1e-8).setMorse(0.15, 30.0, 2.7) },
574 { iListInput().setFene(0.4, 5.0) }
577 //! Constants for Quartic Angles
578 const real cQuarticAngles[5] = { 1.1, 2.3, 4.6, 7.8, 9.2 };
580 //! Function types for testing angles. Add new terms at the end.
581 std::vector<iListInput> c_InputAngles = {
582 { iListInput(2e-3, 1e-8).setHarmonic(F_ANGLES, 100.0, 50.0) },
583 { iListInput(2e-3, 1e-8).setHarmonic(F_ANGLES, 100.15, 50.0, 95.0, 30.0) },
584 { iListInput(8e-3, 1e-8).setHarmonic(F_G96ANGLES, 100.0, 50.0) },
585 { iListInput(8e-3, 1e-8).setHarmonic(F_G96ANGLES, 100.0, 50.0, 95.0, 30.0) },
586 { iListInput().setLinearAngle(50.0, 0.4) },
587 { iListInput().setLinearAngle(50.0, 0.4, 40.0, 0.6) },
588 { iListInput(2e-6, 1e-8).setCrossBondBonds(0.8, 0.7, 45.0) },
589 { iListInput(3e-6, 1e-8).setCrossBondAngles(0.8, 0.7, 0.3, 45.0) },
590 { iListInput(2e-2, 1e-8).setUreyBradley(950.0, 46.0, 0.3, 5.0) },
591 { iListInput(2e-2, 1e-8).setUreyBradley(100.0, 45.0, 0.3, 5.0, 90.0, 47.0, 0.32, 7.0) },
592 { iListInput(2e-3, 1e-8).setQuarticAngles(87.0, cQuarticAngles) }
595 //! Constants for Ryckaert-Bellemans A
596 const real rbcA[NR_RBDIHS] = { -5.35, 13.6, 8.4, -16.7, 0.3, 12.4 };
598 //! Constants for Ryckaert-Bellemans B
599 const real rbcB[NR_RBDIHS] = { -6.35, 12.6, 8.1, -10.7, 0.9, 15.4 };
601 //! Constants for Ryckaert-Bellemans without FEP
602 const real rbc[NR_RBDIHS] = { -7.35, 13.6, 8.4, -16.7, 1.3, 12.4 };
604 //! Function types for testing dihedrals. Add new terms at the end.
605 std::vector<iListInput> c_InputDihs = {
606 { iListInput(5e-4, 1e-8).setPDihedrals(F_PDIHS, -100.0, 10.0, 2, -80.0, 20.0) },
607 { iListInput(1e-4, 1e-8).setPDihedrals(F_PDIHS, -105.0, 15.0, 2) },
608 { iListInput(2e-4, 1e-8).setHarmonic(F_IDIHS, 100.0, 50.0) },
609 { iListInput(2e-4, 1e-8).setHarmonic(F_IDIHS, 100.15, 50.0, 95.0, 30.0) },
610 { iListInput(4e-4, 1e-8).setRbDihedrals(rbcA, rbcB) },
611 { iListInput(4e-4, 1e-8).setRbDihedrals(rbc) }
614 //! Function types for testing polarization. Add new terms at the end.
615 std::vector<iListInput> c_InputPols = {
616 { iListInput(2e-5, 1e-8).setPolarization(0.12) },
617 { iListInput(2e-3, 1e-8).setAnharmPolarization(0.0013, 0.02, 1235.6) },
618 { iListInput(1.4e-3, 1e-8).setTholePolarization(0.26, 0.07, 0.09, 1.6) },
619 { iListInput(2e-3, 1e-8).setWaterPolarization(0.001, 0.0012, 0.0016, 0.095, 0.15, 0.02) },
622 //! Function types for testing polarization. Add new terms at the end.
623 std::vector<iListInput> c_InputRestraints = {
624 { iListInput(1e-4, 1e-8).setPDihedrals(F_ANGRES, -100.0, 10.0, 2, -80.0, 20.0) },
625 { iListInput(1e-4, 1e-8).setPDihedrals(F_ANGRES, -105.0, 15.0, 2) },
626 { iListInput(1e-4, 1e-8).setPDihedrals(F_ANGRESZ, -100.0, 10.0, 2, -80.0, 20.0) },
627 { iListInput(1e-4, 1e-8).setPDihedrals(F_ANGRESZ, -105.0, 15.0, 2) },
628 { iListInput(2e-3, 1e-8).setHarmonic(F_RESTRANGLES, 100.0, 50.0) },
629 { iListInput(2e-3, 1e-8).setHarmonic(F_RESTRANGLES, 100.0, 50.0, 110.0, 45.0) }
632 //! Coordinates for testing
633 std::vector<std::vector<gmx::RVec>> c_coordinatesForTests = {
634 { { 0.0, 0.0, 0.0 }, { 0.0, 0.0, 0.2 }, { 0.005, 0.0, 0.1 }, { -0.001, 0.1, 0.0 } },
635 { { 0.5, 0.0, 0.0 }, { 0.5, 0.0, 0.15 }, { 0.5, 0.07, 0.22 }, { 0.5, 0.18, 0.22 } },
636 { { -0.1143, -0.0282, 0.0 }, { 0.0, 0.0434, 0.0 }, { 0.1185, -0.0138, 0.0 }, { -0.0195, 0.1498, 0.0 } }
639 //! PBC values for testing
640 std::vector<int> c_pbcForTests = { epbcNONE, epbcXY, epbcXYZ };
642 // Those tests give errors with the intel compiler and nothing else, so we disable them only there.
643 #ifndef __INTEL_COMPILER
644 INSTANTIATE_TEST_CASE_P(Bond,
646 ::testing::Combine(::testing::ValuesIn(c_InputBonds),
647 ::testing::ValuesIn(c_coordinatesForTests),
648 ::testing::ValuesIn(c_pbcForTests)));
650 INSTANTIATE_TEST_CASE_P(Angle,
652 ::testing::Combine(::testing::ValuesIn(c_InputAngles),
653 ::testing::ValuesIn(c_coordinatesForTests),
654 ::testing::ValuesIn(c_pbcForTests)));
656 INSTANTIATE_TEST_CASE_P(Dihedral,
658 ::testing::Combine(::testing::ValuesIn(c_InputDihs),
659 ::testing::ValuesIn(c_coordinatesForTests),
660 ::testing::ValuesIn(c_pbcForTests)));
662 INSTANTIATE_TEST_CASE_P(Polarize,
664 ::testing::Combine(::testing::ValuesIn(c_InputPols),
665 ::testing::ValuesIn(c_coordinatesForTests),
666 ::testing::ValuesIn(c_pbcForTests)));
668 INSTANTIATE_TEST_CASE_P(Restraints,
670 ::testing::Combine(::testing::ValuesIn(c_InputRestraints),
671 ::testing::ValuesIn(c_coordinatesForTests),
672 ::testing::ValuesIn(c_pbcForTests)));