1 /*
2 * This file is part of the GROMACS molecular simulation package.
3 *
4 * Copyright (c) 2018,2019,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 Defines the runner for CUDA version of SETTLE.
37 *
38 * \author Artem Zhmurov <zhmurov@gmail.com>
39 * \ingroup module_mdlib
40 */
41 #include "gmxpre.h"
42
43 #include "settletestrunners.h"
44
45 #include "config.h"
46
47 #include <assert.h>
48
49 #include <cmath>
50
51 #include <algorithm>
52 #include <vector>
53
54 #include "gromacs/gpu_utils/devicebuffer.cuh"
55 #include "gromacs/hardware/device_information.h"
56 #include "gromacs/mdlib/settle_gpu.cuh"
57 #include "gromacs/utility/unique_cptr.h"
58
59 #include "testutils/test_device.h"
60
61 namespace gmx
62 {
63 namespace test
64 {
65
applySettle(SettleTestData * testData,const t_pbc pbc,const bool updateVelocities,const bool calcVirial,const std::string &)66 void SettleDeviceTestRunner::applySettle(SettleTestData* testData,
67 const t_pbc pbc,
68 const bool updateVelocities,
69 const bool calcVirial,
70 const std::string& /* testDescription */)
71 {
72 // These should never fail since this function should only be called if CUDA is enabled and
73 // there is a CUDA-capable device available.
74 GMX_RELEASE_ASSERT(GMX_GPU_CUDA, "CUDA version of SETTLE was called from non-CUDA build.");
75
76 const DeviceContext& deviceContext = testDevice_.deviceContext();
77 const DeviceStream& deviceStream = testDevice_.deviceStream();
78 setActiveDevice(testDevice_.deviceInfo());
79
80 auto settleGpu = std::make_unique<SettleGpu>(testData->mtop_, deviceContext, deviceStream);
81
82 settleGpu->set(*testData->idef_);
83 PbcAiuc pbcAiuc;
84 setPbcAiuc(pbc.ndim_ePBC, pbc.box, &pbcAiuc);
85
86 int numAtoms = testData->numAtoms_;
87
88 float3 *d_x, *d_xp, *d_v;
89
90 float3* h_x = (float3*)(as_rvec_array(testData->x_.data()));
91 float3* h_xp = (float3*)(as_rvec_array(testData->xPrime_.data()));
92 float3* h_v = (float3*)(as_rvec_array(testData->v_.data()));
93
94 allocateDeviceBuffer(&d_x, numAtoms, deviceContext);
95 allocateDeviceBuffer(&d_xp, numAtoms, deviceContext);
96 allocateDeviceBuffer(&d_v, numAtoms, deviceContext);
97
98 copyToDeviceBuffer(&d_x, (float3*)h_x, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
99 copyToDeviceBuffer(&d_xp, (float3*)h_xp, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
100 if (updateVelocities)
101 {
102 copyToDeviceBuffer(&d_v, (float3*)h_v, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
103 }
104 settleGpu->apply(d_x, d_xp, updateVelocities, d_v, testData->reciprocalTimeStep_, calcVirial,
105 testData->virial_, pbcAiuc);
106
107 copyFromDeviceBuffer((float3*)h_xp, &d_xp, 0, numAtoms, deviceStream, GpuApiCallBehavior::Sync, nullptr);
108 if (updateVelocities)
109 {
110 copyFromDeviceBuffer((float3*)h_v, &d_v, 0, numAtoms, deviceStream,
111 GpuApiCallBehavior::Sync, nullptr);
112 }
113
114 freeDeviceBuffer(&d_x);
115 freeDeviceBuffer(&d_xp);
116 freeDeviceBuffer(&d_v);
117 }
118
119 } // namespace test
120 } // namespace gmx
121