816a545dfb3067ecc7c9db643aed609818571c02
[alexxy/gromacs.git] / api / nblib / tests / nbkernelsystem.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
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 topology setup 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  * \author Artem Zhmurov <zhmurov@gmail.com>
44  */
45 #include <gtest/gtest.h>
46
47 #include "gromacs/topology/exclusionblocks.h"
48 #include "nblib/forcecalculator.h"
49 #include "nblib/gmxsetup.h"
50 #include "nblib/integrator.h"
51 #include "nblib/tests/testhelpers.h"
52 #include "nblib/tests/testsystems.h"
53 #include "nblib/topology.h"
54
55 namespace nblib
56 {
57 namespace test
58 {
59 namespace
60 {
61
62 // This is defined in src/gromacs/mdtypes/forcerec.h but there is also a
63 // legacy C6 macro defined there that conflicts with the nblib C6 type.
64 // Todo: Once that C6 has been refactored into a regular function, this
65 //       file can just include forcerec.h
66 //! Macro to set Van der Waals interactions to atoms
67 #define SET_CGINFO_HAS_VDW(cgi) (cgi) = ((cgi) | (1 << 23))
68
69 TEST(NBlibTest, SpcMethanolForcesAreCorrect)
70 {
71     auto options        = NBKernelOptions();
72     options.nbnxmSimd   = SimdKernels::SimdNo;
73     options.coulombType = CoulombType::Cutoff;
74
75     SpcMethanolSimulationStateBuilder spcMethanolSystemBuilder;
76
77     auto simState        = spcMethanolSystemBuilder.setupSimulationState();
78     auto forceCalculator = ForceCalculator(simState, options);
79
80     gmx::ArrayRef<Vec3> forces(simState.forces());
81     ASSERT_NO_THROW(forceCalculator.compute(simState.coordinates(), forces));
82
83     /* Use higher-than-usual tolerance for forces. Some of the particles in the test systems are
84      * very close to each other, and, for example, the distance between the first two particles
85      * is approx. 0.13 and already has relative uncertainty around 1e-6. */
86     gmx::test::FloatingPointTolerance forceTolerance(1.0e-5, 1.0e-9, 1e-4, 1.0e-9, 1000, 1000, true);
87
88     Vector3DTest forcesOutputTest(forceTolerance);
89     forcesOutputTest.testVectors(forces, "SPC-methanol forces");
90 }
91
92 TEST(NBlibTest, ExpectedNumberOfForces)
93 {
94     auto options      = NBKernelOptions();
95     options.nbnxmSimd = SimdKernels::SimdNo;
96
97     SpcMethanolSimulationStateBuilder spcMethanolSystemBuilder;
98
99     auto simState        = spcMethanolSystemBuilder.setupSimulationState();
100     auto forceCalculator = ForceCalculator(simState, options);
101
102     gmx::ArrayRef<Vec3> forces(simState.forces());
103     forceCalculator.compute(simState.coordinates(), forces);
104     EXPECT_EQ(simState.topology().numParticles(), forces.size());
105 }
106
107 TEST(NBlibTest, CanIntegrateSystem)
108 {
109     auto options          = NBKernelOptions();
110     options.nbnxmSimd     = SimdKernels::SimdNo;
111     options.numIterations = 1;
112
113     SpcMethanolSimulationStateBuilder spcMethanolSystemBuilder;
114
115     auto simState        = spcMethanolSystemBuilder.setupSimulationState();
116     auto forceCalculator = ForceCalculator(simState, options);
117
118     LeapFrog integrator(simState.topology(), simState.box());
119
120     for (int iter = 0; iter < options.numIterations; iter++)
121     {
122         gmx::ArrayRef<Vec3> forces(simState.forces());
123         forceCalculator.compute(simState.coordinates(), simState.forces());
124         EXPECT_NO_THROW(integrator.integrate(1.0, simState.coordinates(), simState.velocities(),
125                                              simState.forces()));
126     }
127 }
128
129 /*!
130  * Check if the following aspects of the ForceCalculator and
131  * LeapFrog (integrator) work as expected:
132  *
133  * 1. Calling the ForceCalculator::compute() function makes no change
134  *    to the internal representation of the system. Calling it repeatedly
135  *    without update should return the same vector of forces.
136  *
137  * 2. Once the LeapFrog objects integrates for the given time using the
138  *    forces, there the coordinates in SimulationState must change.
139  *    Calling the compute() function must now generate a new set of forces.
140  *
141  */
142 TEST(NBlibTest, UpdateChangesForces)
143 {
144     auto options          = NBKernelOptions();
145     options.nbnxmSimd     = SimdKernels::SimdNo;
146     options.numIterations = 1;
147
148     SpcMethanolSimulationStateBuilder spcMethanolSystemBuilder;
149
150     auto simState        = spcMethanolSystemBuilder.setupSimulationState();
151     auto forceCalculator = ForceCalculator(simState, options);
152
153     LeapFrog integrator(simState.topology(), simState.box());
154
155     // step 1
156     gmx::ArrayRef<Vec3> forces(simState.forces());
157     forceCalculator.compute(simState.coordinates(), simState.forces());
158
159     // copy computed forces to another array
160     std::vector<Vec3> forces_1(forces.size());
161     std::copy(forces.begin(), forces.end(), begin(forces_1));
162
163     // zero original force buffer
164     zeroCartesianArray(forces);
165
166     // check if forces change without update step
167     forceCalculator.compute(simState.coordinates(), forces);
168
169     // check if forces change without update
170     for (size_t i = 0; i < forces_1.size(); i++)
171     {
172         for (int j = 0; j < dimSize; j++)
173         {
174             EXPECT_EQ(forces[i][j], forces_1[i][j]);
175         }
176     }
177
178     // update
179     integrator.integrate(1.0, simState.coordinates(), simState.velocities(), simState.forces());
180
181     // zero original force buffer
182     zeroCartesianArray(forces);
183
184     // step 2
185     forceCalculator.compute(simState.coordinates(), forces);
186     std::vector<Vec3> forces_2(forces.size());
187     std::copy(forces.begin(), forces.end(), begin(forces_2));
188
189     // check if forces change after update
190     for (size_t i = 0; i < forces_1.size(); i++)
191     {
192         for (int j = 0; j < dimSize; j++)
193         {
194             EXPECT_NE(forces_1[i][j], forces_2[i][j]);
195         }
196     }
197 }
198
199 TEST(NBlibTest, ArgonForcesAreCorrect)
200 {
201     auto options        = NBKernelOptions();
202     options.nbnxmSimd   = SimdKernels::SimdNo;
203     options.coulombType = CoulombType::Cutoff;
204
205     ArgonSimulationStateBuilder argonSystemBuilder;
206
207     auto simState        = argonSystemBuilder.setupSimulationState();
208     auto forceCalculator = ForceCalculator(simState, options);
209
210     gmx::ArrayRef<Vec3> testForces(simState.forces());
211     forceCalculator.compute(simState.coordinates(), simState.forces());
212
213     Vector3DTest forcesOutputTest;
214     forcesOutputTest.testVectors(testForces, "Argon forces");
215 }
216
217 } // namespace
218 } // namespace test
219 } // namespace nblib