Set temperature scaling groups upon construction of integrator
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / leapfrogtestrunners_gpu.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,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 Runner for GPU version of the integrator
37  *
38  * Handles GPU data management and actual numerical integration.
39  *
40  * \author Artem Zhmurov <zhmurov@gmail.com>
41  * \ingroup module_mdlib
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #include <gtest/gtest.h>
48
49 #include "leapfrogtestrunners.h"
50
51 #if GMX_GPU_CUDA
52 #    include "gromacs/gpu_utils/devicebuffer.cuh"
53 #endif
54 #if GMX_GPU_SYCL
55 #    include "gromacs/gpu_utils/devicebuffer_sycl.h"
56 #endif
57
58 #if HAVE_GPU_LEAPFROG
59 #    include "gromacs/mdlib/leapfrog_gpu.h"
60 #endif
61
62 #include "gromacs/hardware/device_information.h"
63 #include "gromacs/mdlib/stat.h"
64
65 namespace gmx
66 {
67 namespace test
68 {
69
70 #if HAVE_GPU_LEAPFROG
71 void LeapFrogDeviceTestRunner::integrate(LeapFrogTestData* testData, int numSteps)
72 {
73     const DeviceContext& deviceContext = testDevice_.deviceContext();
74     const DeviceStream&  deviceStream  = testDevice_.deviceStream();
75     setActiveDevice(testDevice_.deviceInfo());
76
77     int numAtoms = testData->numAtoms_;
78
79     static_assert(sizeof(float3) == sizeof(*testData->x_.data()), "Incompatible types");
80
81     float3* h_x  = reinterpret_cast<float3*>(testData->x_.data());
82     float3* h_xp = reinterpret_cast<float3*>(testData->xPrime_.data());
83     float3* h_v  = reinterpret_cast<float3*>(testData->v_.data());
84     float3* h_f  = reinterpret_cast<float3*>(testData->f_.data());
85
86     DeviceBuffer<float3> d_x, d_xp, d_v, d_f;
87
88     allocateDeviceBuffer(&d_x, numAtoms, deviceContext);
89     allocateDeviceBuffer(&d_xp, numAtoms, deviceContext);
90     allocateDeviceBuffer(&d_v, numAtoms, deviceContext);
91     allocateDeviceBuffer(&d_f, numAtoms, deviceContext);
92
93     copyToDeviceBuffer(&d_x, h_x, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
94     copyToDeviceBuffer(&d_xp, h_xp, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
95     copyToDeviceBuffer(&d_v, h_v, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
96     copyToDeviceBuffer(&d_f, h_f, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
97
98     auto integrator =
99             std::make_unique<LeapFrogGpu>(deviceContext, deviceStream, testData->numTCoupleGroups_);
100
101     integrator->set(testData->numAtoms_, testData->inverseMasses_.data(), testData->mdAtoms_.cTC);
102
103     bool doTempCouple = testData->numTCoupleGroups_ > 0;
104     for (int step = 0; step < numSteps; step++)
105     {
106         // This follows the logic of the CPU-based implementation
107         bool doPressureCouple = testData->doPressureCouple_
108                                 && do_per_step(step + testData->inputRecord_.nstpcouple - 1,
109                                                testData->inputRecord_.nstpcouple);
110         integrator->integrate(d_x,
111                               d_xp,
112                               d_v,
113                               d_f,
114                               testData->timestep_,
115                               doTempCouple,
116                               testData->kineticEnergyData_.tcstat,
117                               doPressureCouple,
118                               testData->dtPressureCouple_,
119                               testData->velocityScalingMatrix_);
120     }
121
122     copyFromDeviceBuffer(h_xp, &d_x, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
123     copyFromDeviceBuffer(h_v, &d_v, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
124
125     freeDeviceBuffer(&d_x);
126     freeDeviceBuffer(&d_xp);
127     freeDeviceBuffer(&d_v);
128     freeDeviceBuffer(&d_f);
129 }
130
131 #else // HAVE_GPU_LEAPFROG
132
133 void LeapFrogDeviceTestRunner::integrate(LeapFrogTestData* /* testData */, int /* numSteps */)
134 {
135     GMX_UNUSED_VALUE(testDevice_);
136     FAIL() << "Dummy Leap-Frog GPU function was called instead of the real one.";
137 }
138
139 #endif // HAVE_GPU_LEAPFROG
140
141 } // namespace test
142 } // namespace gmx