Introduce nblib nonbonded force calculator impl class
[alexxy/gromacs.git] / api / nblib / tests / gmxcalculator.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 utility 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 <gtest/gtest.h>
45
46 #include "nblib/gmxcalculatorcpu.h"
47 #include "nblib/kerneloptions.h"
48 #include "nblib/simulationstate.h"
49 #include "nblib/tests/testhelpers.h"
50 #include "nblib/tests/testsystems.h"
51 #include "gromacs/utility/arrayref.h"
52
53 namespace nblib
54 {
55 namespace test
56 {
57 namespace
58 {
59 TEST(NBlibTest, GmxForceCalculatorCanCompute)
60 {
61     ArgonSimulationStateBuilder argonSystemBuilder(fftypes::GROMOS43A1);
62     SimulationState             simState = argonSystemBuilder.setupSimulationState();
63     NBKernelOptions             options  = NBKernelOptions();
64     options.nbnxmSimd                    = SimdKernels::SimdNo;
65
66     std::unique_ptr<GmxNBForceCalculatorCpu> gmxForceCalculator =
67             setupGmxForceCalculatorCpu(simState.topology(), options);
68     gmxForceCalculator->updatePairlist(simState.coordinates(), simState.box());
69
70     EXPECT_NO_THROW(gmxForceCalculator->compute(simState.coordinates(), simState.box(), simState.forces()));
71 }
72
73 TEST(NBlibTest, ArgonVirialsAreCorrect)
74 {
75     ArgonSimulationStateBuilder argonSystemBuilder(fftypes::OPLSA);
76     SimulationState             simState = argonSystemBuilder.setupSimulationState();
77     NBKernelOptions             options  = NBKernelOptions();
78     options.nbnxmSimd                    = SimdKernels::SimdNo;
79     std::unique_ptr<GmxNBForceCalculatorCpu> gmxForceCalculator =
80             setupGmxForceCalculatorCpu(simState.topology(), options);
81     gmxForceCalculator->updatePairlist(simState.coordinates(), simState.box());
82
83     std::vector<real> virialArray(9, 0.0);
84
85     gmxForceCalculator->compute(simState.coordinates(), simState.box(), simState.forces(), virialArray);
86
87     RefDataChecker virialsOutputTest(1e-7);
88     virialsOutputTest.testArrays<real>(virialArray, "Virials");
89 }
90
91 TEST(NBlibTest, ArgonEnergiesAreCorrect)
92 {
93     ArgonSimulationStateBuilder argonSystemBuilder(fftypes::OPLSA);
94     SimulationState             simState = argonSystemBuilder.setupSimulationState();
95     NBKernelOptions             options  = NBKernelOptions();
96     options.nbnxmSimd                    = SimdKernels::SimdNo;
97     std::unique_ptr<GmxNBForceCalculatorCpu> gmxForceCalculator =
98             setupGmxForceCalculatorCpu(simState.topology(), options);
99     gmxForceCalculator->updatePairlist(simState.coordinates(), simState.box());
100
101     // number of energy kinds is 5: COULSR, LJSR, BHAMSR, COUL14, LJ14,
102     std::vector<real> energies(5, 0.0);
103
104     gmxForceCalculator->compute(
105             simState.coordinates(), simState.box(), simState.forces(), gmx::ArrayRef<real>{}, energies);
106
107     RefDataChecker energiesOutputTest(5e-5);
108     energiesOutputTest.testArrays<real>(energies, "Argon energies");
109 }
110
111 TEST(NBlibTest, SpcMethanolEnergiesAreCorrect)
112 {
113     SpcMethanolSimulationStateBuilder spcMethanolSystemBuilder;
114     SimulationState                   simState = spcMethanolSystemBuilder.setupSimulationState();
115     NBKernelOptions                   options  = NBKernelOptions();
116     options.nbnxmSimd                          = SimdKernels::SimdNo;
117     std::unique_ptr<GmxNBForceCalculatorCpu> gmxForceCalculator =
118             setupGmxForceCalculatorCpu(simState.topology(), options);
119     gmxForceCalculator->updatePairlist(simState.coordinates(), simState.box());
120
121     // number of energy kinds is 5: COULSR, LJSR, BHAMSR, COUL14, LJ14,
122     std::vector<real> energies(5, 0.0);
123
124     gmxForceCalculator->compute(
125             simState.coordinates(), simState.box(), simState.forces(), gmx::ArrayRef<real>{}, energies);
126
127     RefDataChecker energiesOutputTest(5e-5);
128     energiesOutputTest.testArrays<real>(energies, "SPC-methanol energies");
129 }
130
131 } // namespace
132 } // namespace test
133 } // namespace nblib