Add several more listed forces interaction types to nblib
[alexxy/gromacs.git] / api / nblib / listed_forces / tests / typetests.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * This implements basic nblib box tests
38  *
39  * \author Victor Holanda <victor.holanda@cscs.ch>
40  * \author Joe Jordan <ejjordan@kth.se>
41  * \author Prashanth Kanduri <kanduri@cscs.ch>
42  * \author Sebastian Keller <keller@cscs.ch>
43  */
44 #include "gromacs/utility/arrayref.h"
45 #include "nblib/listed_forces/dataflow.hpp"
46 #include "nblib/listed_forces/tests/listedtesthelpers.h"
47 #include "nblib/tests/testhelpers.h"
48
49 #include "testutils/refdata.h"
50 #include "testutils/testasserts.h"
51
52 namespace nblib
53 {
54
55 //! Coordinates for testing
56 static const std::vector<gmx::RVec> c_coordinatesForDihTests = { { 0.0, 0.0, 0.0 },
57                                                                  { 0.0, 0.0, 0.2 },
58                                                                  { 0.005, 0.0, 0.1 },
59                                                                  { -0.001, 0.1, 0.0 } };
60
61 //! Coordinates for testing angles
62 static const std::vector<gmx::RVec> c_coordinatesForAngleTests = { { 1.382, 1.573, 1.482 },
63                                                                    { 1.281, 1.559, 1.596 },
64                                                                    { 1.292, 1.422, 1.663 } };
65
66 //! Coordinates for testing bonds
67 static const std::vector<gmx::RVec> c_coordinatesForBondTests = { { 1.382, 1.573, 1.482 },
68                                                                   { 1.281, 1.559, 1.596 } };
69
70 //! Function types for testing Harmonic bonds
71 static const std::vector<HarmonicBondType> c_InputHarmonicBonds = { { HarmonicBondType(500.0, 0.15) } };
72
73 //! Function types for testing G96 bonds
74 static const std::vector<G96BondType> c_InputG96Bonds = { { G96BondType(50.0, 0.15) } };
75
76 //! Function types for testing cubic bonds
77 static const std::vector<CubicBondType> c_InputCubicBonds = { { CubicBondType(50.0, 2.0, 0.16) } };
78
79 //! Function types for testing Morse bonds
80 static const std::vector<MorseBondType> c_InputMorseBonds = { { MorseBondType(30.0, 2.7, 0.15) } };
81
82 //! Function types for testing FENE bonds
83 static const std::vector<FENEBondType> c_InputFeneBonds = { { FENEBondType(5.0, 0.4) } };
84
85 //! Function types for testing Harmonic angles
86 static const std::vector<HarmonicAngle> c_InputHarmonicAngles = { { HarmonicAngle(50.0, Degrees(100)) } };
87
88 //! Function types for testing Linear angles
89 static const std::vector<LinearAngle> c_InputLinearAngles = { { LinearAngle(50.0, 0.4) } };
90
91 //! Function types for testing G96 angles
92 static const std::vector<G96Angle> c_InputG96Angles = { { G96Angle(50.0, Degrees(100)) } };
93
94 //! Function types for testing Restricted angles
95 static const std::vector<RestrictedAngle> c_InputRestrictedAngles = { { RestrictedAngle(50.0, Degrees(100)) } };
96
97 //! Function types for testing Quartic angles
98 static const std::vector<QuarticAngle> c_InputQuarticAngles = {
99     { QuarticAngle(1.1, 2.3, 4.6, 7.8, 9.2, Degrees(87)) }
100 };
101
102 //! Function types for testing cross bond-bond interaction
103 static const std::vector<CrossBondBond> c_InputCrossBondBond = { { CrossBondBond(45.0, 0.8, 0.7) } };
104
105 //! Function types for testing cross bond-angle interaction
106 static const std::vector<CrossBondAngle> c_InputCrossBondAngle = { { CrossBondAngle(45.0, 0.8, 0.7, 0.3) } };
107
108 //! Function types for testing dihedrals
109 static const std::vector<ProperDihedral> c_InputDihs = { { ProperDihedral(Degrees(-105.0), 15.0, 2) } };
110
111 template<class Interaction, std::enable_if_t<Contains<Interaction, SupportedListedTypes>{}>* = nullptr>
112 void checkForcesAndEnergiesWithRefData(std::vector<Interaction> input, gmx::ArrayRef<const gmx::RVec> x)
113 {
114     auto                   indices = indexVector<Interaction>();
115     PbcHolder              pbcHolder(PbcType::Xyz, Box(1.5));
116     test::RefDataChecker   refDataChecker(1e-4);
117     std::vector<gmx::RVec> forces(x.size(), gmx::RVec{ 0, 0, 0 });
118
119     auto energy = computeForces(gmx::ArrayRef<const InteractionIndex<Interaction>>(indices),
120                                 gmx::ArrayRef<const Interaction>(input),
121                                 x,
122                                 &forces,
123                                 pbcHolder);
124
125     refDataChecker.testReal(energy, "Epot");
126     refDataChecker.testArrays<gmx::RVec>(forces, "forces");
127 }
128
129 TEST(FourCenter, ListedForcesProperDihedralTest)
130 {
131     checkForcesAndEnergiesWithRefData(c_InputDihs, c_coordinatesForDihTests);
132 }
133
134 TEST(ThreeCenter, ListedForcesG96AngleTest)
135 {
136     checkForcesAndEnergiesWithRefData(c_InputG96Angles, c_coordinatesForAngleTests);
137 }
138
139 TEST(ThreeCenter, ListedForcesHarmonicAngleTest)
140 {
141     checkForcesAndEnergiesWithRefData(c_InputHarmonicAngles, c_coordinatesForAngleTests);
142 }
143
144 TEST(ThreeCenter, ListedForcesLinearAngleTest)
145 {
146     checkForcesAndEnergiesWithRefData(c_InputLinearAngles, c_coordinatesForAngleTests);
147 }
148
149 TEST(ThreeCenter, ListedForcesCrossBondBondTest)
150 {
151     checkForcesAndEnergiesWithRefData(c_InputCrossBondBond, c_coordinatesForAngleTests);
152 }
153
154 TEST(ThreeCenter, ListedForcesCrossBondAngleTest)
155 {
156     checkForcesAndEnergiesWithRefData(c_InputCrossBondAngle, c_coordinatesForAngleTests);
157 }
158
159 TEST(ThreeCenter, ListedForcesQuarticAngleTest)
160 {
161     checkForcesAndEnergiesWithRefData(c_InputQuarticAngles, c_coordinatesForAngleTests);
162 }
163
164 TEST(ThreeCenter, ListedForcesRestrictedAngleTest)
165 {
166     checkForcesAndEnergiesWithRefData(c_InputRestrictedAngles, c_coordinatesForAngleTests);
167 }
168
169 TEST(TwoCenter, ListedForcesHarmonicBondTest)
170 {
171     checkForcesAndEnergiesWithRefData(c_InputHarmonicBonds, c_coordinatesForBondTests);
172 }
173
174 TEST(TwoCenter, ListedForcesG96BondTest)
175 {
176     checkForcesAndEnergiesWithRefData(c_InputG96Bonds, c_coordinatesForBondTests);
177 }
178
179 TEST(TwoCenter, ListedForcesCubicBondTest)
180 {
181     checkForcesAndEnergiesWithRefData(c_InputCubicBonds, c_coordinatesForBondTests);
182 }
183
184 TEST(TwoCenter, ListedForcesMorseBondTest)
185 {
186     checkForcesAndEnergiesWithRefData(c_InputMorseBonds, c_coordinatesForBondTests);
187 }
188
189 TEST(TwoCenter, ListedForcesFeneBondTest)
190 {
191     checkForcesAndEnergiesWithRefData(c_InputFeneBonds, c_coordinatesForBondTests);
192 }
193
194 } // namespace nblib