Unify nbnxn_gpu_init_x_to_nbat_x
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_cuda_data_mgmt.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \file
37  *  \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
38  *
39  *  \author Szilard Pall <pall.szilard@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include <assert.h>
44 #include <stdarg.h>
45 #include <stdio.h>
46 #include <stdlib.h>
47
48 // TODO We would like to move this down, but the way NbnxmGpu
49 //      is currently declared means this has to be before gpu_types.h
50 #include "nbnxm_cuda_types.h"
51
52 // TODO Remove this comment when the above order issue is resolved
53 #include "gromacs/gpu_utils/cudautils.cuh"
54 #include "gromacs/gpu_utils/device_context.h"
55 #include "gromacs/gpu_utils/gpu_utils.h"
56 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
57 #include "gromacs/gpu_utils/pmalloc.h"
58 #include "gromacs/hardware/device_information.h"
59 #include "gromacs/hardware/device_management.h"
60 #include "gromacs/math/vectypes.h"
61 #include "gromacs/mdlib/force_flags.h"
62 #include "gromacs/mdtypes/interaction_const.h"
63 #include "gromacs/mdtypes/md_enums.h"
64 #include "gromacs/nbnxm/atomdata.h"
65 #include "gromacs/nbnxm/gpu_data_mgmt.h"
66 #include "gromacs/nbnxm/gridset.h"
67 #include "gromacs/nbnxm/nbnxm.h"
68 #include "gromacs/nbnxm/nbnxm_gpu.h"
69 #include "gromacs/nbnxm/nbnxm_gpu_data_mgmt.h"
70 #include "gromacs/nbnxm/pairlistsets.h"
71 #include "gromacs/pbcutil/ishift.h"
72 #include "gromacs/timing/gpu_timing.h"
73 #include "gromacs/utility/basedefinitions.h"
74 #include "gromacs/utility/cstringutil.h"
75 #include "gromacs/utility/fatalerror.h"
76 #include "gromacs/utility/real.h"
77 #include "gromacs/utility/smalloc.h"
78
79 #include "nbnxm_cuda.h"
80
81 namespace Nbnxm
82 {
83
84 /* This is a heuristically determined parameter for the Kepler
85  * and Maxwell architectures for the minimum size of ci lists by multiplying
86  * this constant with the # of multiprocessors on the current device.
87  * Since the maximum number of blocks per multiprocessor is 16, the ideal
88  * count for small systems is 32 or 48 blocks per multiprocessor. Because
89  * there is a bit of fluctuations in the generated block counts, we use
90  * a target of 44 instead of the ideal value of 48.
91  */
92 static unsigned int gpu_min_ci_balanced_factor = 44;
93
94 void gpu_init_platform_specific(NbnxmGpu* /* nb */)
95 {
96     /* set the kernel type for the current GPU */
97     /* pick L1 cache configuration */
98     cuda_set_cacheconfig();
99 }
100
101 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
102 {
103     NBAtomData*         adat        = nb->atdat;
104     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
105
106     /* only if we have a dynamic box */
107     if (nbatom->bDynamicBox || !adat->shiftVecUploaded)
108     {
109         static_assert(sizeof(adat->shiftVec[0]) == sizeof(nbatom->shift_vec[0]),
110                       "Sizes of host- and device-side shift vectors should be the same.");
111         copyToDeviceBuffer(&adat->shiftVec,
112                            reinterpret_cast<const Float3*>(nbatom->shift_vec.data()),
113                            0,
114                            SHIFTS,
115                            localStream,
116                            GpuApiCallBehavior::Async,
117                            nullptr);
118         adat->shiftVecUploaded = true;
119     }
120 }
121
122 void gpu_free(NbnxmGpu* nb)
123 {
124     if (nb == nullptr)
125     {
126         return;
127     }
128
129     delete nb->timers;
130     sfree(nb->timings);
131
132     NBAtomData* atdat   = nb->atdat;
133     NBParamGpu* nbparam = nb->nbparam;
134
135     if (nbparam->elecType == ElecType::EwaldTab || nbparam->elecType == ElecType::EwaldTabTwin)
136     {
137         destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
138     }
139
140     if (!useLjCombRule(nb->nbparam->vdwType))
141     {
142         destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
143     }
144
145     if (nbparam->vdwType == VdwType::EwaldGeom || nbparam->vdwType == VdwType::EwaldLB)
146     {
147         destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
148     }
149
150     freeDeviceBuffer(&atdat->shiftVec);
151     freeDeviceBuffer(&atdat->fShift);
152
153     freeDeviceBuffer(&atdat->eLJ);
154     freeDeviceBuffer(&atdat->eElec);
155
156     freeDeviceBuffer(&atdat->f);
157     freeDeviceBuffer(&atdat->xq);
158     if (useLjCombRule(nb->nbparam->vdwType))
159     {
160         freeDeviceBuffer(&atdat->ljComb);
161     }
162     else
163     {
164         freeDeviceBuffer(&atdat->atomTypes);
165     }
166
167     /* Free plist */
168     auto* plist = nb->plist[InteractionLocality::Local];
169     freeDeviceBuffer(&plist->sci);
170     freeDeviceBuffer(&plist->cj4);
171     freeDeviceBuffer(&plist->imask);
172     freeDeviceBuffer(&plist->excl);
173     delete plist;
174     if (nb->bUseTwoStreams)
175     {
176         auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
177         freeDeviceBuffer(&plist_nl->sci);
178         freeDeviceBuffer(&plist_nl->cj4);
179         freeDeviceBuffer(&plist_nl->imask);
180         freeDeviceBuffer(&plist_nl->excl);
181         delete plist_nl;
182     }
183
184     /* Free nbst */
185     pfree(nb->nbst.eLJ);
186     nb->nbst.eLJ = nullptr;
187
188     pfree(nb->nbst.eElec);
189     nb->nbst.eElec = nullptr;
190
191     pfree(nb->nbst.fShift);
192     nb->nbst.fShift = nullptr;
193
194     delete atdat;
195     delete nbparam;
196     delete nb;
197
198     if (debug)
199     {
200         fprintf(debug, "Cleaned up CUDA data structures.\n");
201     }
202 }
203
204 int gpu_min_ci_balanced(NbnxmGpu* nb)
205 {
206     return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
207                          : 0;
208 }
209
210 void* gpu_get_xq(NbnxmGpu* nb)
211 {
212     assert(nb);
213
214     return static_cast<void*>(nb->atdat->xq);
215 }
216
217 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
218 {
219     assert(nb);
220
221     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
222 }
223
224 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
225 {
226     assert(nb);
227
228     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fShift);
229 }
230
231 } // namespace Nbnxm