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