xref: /libCEED/backends/cuda-shared/ceed-cuda-shared-basis.c (revision 8ad9a6ff7220956001918baa07d5a82989f396b2)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include <ceed.h>
9 #include <ceed/backend.h>
10 #include <ceed/jit-tools.h>
11 #include <cuda.h>
12 #include <cuda_runtime.h>
13 #include <stdbool.h>
14 #include <stddef.h>
15 
16 #include "../cuda/ceed-cuda-common.h"
17 #include "../cuda/ceed-cuda-compile.h"
18 #include "ceed-cuda-shared.h"
19 
20 //------------------------------------------------------------------------------
21 // Device initalization
22 //------------------------------------------------------------------------------
23 int CeedInit_CudaInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B);
24 int CeedInit_CudaGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
25 int CeedInit_CudaCollocatedGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr);
26 
27 //------------------------------------------------------------------------------
28 // Apply basis
29 //------------------------------------------------------------------------------
30 int CeedBasisApplyTensor_Cuda_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u,
31                                      CeedVector v) {
32   Ceed                   ceed;
33   Ceed_Cuda             *ceed_Cuda;
34   CeedInt                dim, num_comp;
35   const CeedScalar      *d_u;
36   CeedScalar            *d_v;
37   CeedBasis_Cuda_shared *data;
38 
39   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
40   CeedCallBackend(CeedGetData(ceed, &ceed_Cuda));
41   CeedCallBackend(CeedBasisGetData(basis, &data));
42   CeedCallBackend(CeedBasisGetDimension(basis, &dim));
43   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
44 
45   // Get read/write access to u, v
46   if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u));
47   else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode");
48   CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v));
49 
50   // Apply basis operation
51   switch (eval_mode) {
52     case CEED_EVAL_INTERP: {
53       CeedInt P_1d, Q_1d;
54 
55       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
56       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
57       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
58 
59       CeedCallBackend(CeedInit_CudaInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B));
60       void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v};
61 
62       if (dim == 1) {
63         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
64                                                                                                  1));  // avoid >512 total threads
65         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
66         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
67 
68         if (t_mode == CEED_TRANSPOSE) {
69           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
70         } else {
71           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args));
72         }
73       } else if (dim == 2) {
74         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
75         // elems_per_block must be at least 1
76         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
77         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
78         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
79 
80         if (t_mode == CEED_TRANSPOSE) {
81           CeedCallBackend(
82               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
83         } else {
84           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
85         }
86       } else if (dim == 3) {
87         CeedInt elems_per_block = 1;
88         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
89         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
90 
91         if (t_mode == CEED_TRANSPOSE) {
92           CeedCallBackend(
93               CeedRunKernelDimShared_Cuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
94         } else {
95           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args));
96         }
97       }
98     } break;
99     case CEED_EVAL_GRAD: {
100       CeedInt P_1d, Q_1d;
101 
102       CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d));
103       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
104       CeedInt thread_1d = CeedIntMax(Q_1d, P_1d);
105 
106       if (data->d_collo_grad_1d) {
107         CeedCallBackend(CeedInit_CudaCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
108       } else {
109         CeedCallBackend(CeedInit_CudaGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G));
110       }
111       void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v};
112       if (dim == 1) {
113         CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d,
114                                                                                                  1));  // avoid >512 total threads
115         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
116         CeedInt shared_mem      = elems_per_block * thread_1d * sizeof(CeedScalar);
117 
118         if (t_mode == CEED_TRANSPOSE) {
119           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
120         } else {
121           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args));
122         }
123       } else if (dim == 2) {
124         const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8};
125         // elems_per_block must be at least 1
126         CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1);
127         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
128         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
129 
130         if (t_mode == CEED_TRANSPOSE) {
131           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
132         } else {
133           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
134         }
135       } else if (dim == 3) {
136         CeedInt elems_per_block = 1;
137         CeedInt grid            = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
138         CeedInt shared_mem      = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar);
139 
140         if (t_mode == CEED_TRANSPOSE) {
141           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
142         } else {
143           CeedCallBackend(CeedRunKernelDimShared_Cuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args));
144         }
145       }
146     } break;
147     case CEED_EVAL_WEIGHT: {
148       CeedInt Q_1d;
149       CeedInt block_size = 32;
150 
151       CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d));
152       void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v};
153       if (dim == 1) {
154         const CeedInt elems_per_block = block_size / Q_1d;
155         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
156 
157         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args));
158       } else if (dim == 2) {
159         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
160         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
161         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
162 
163         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
164       } else if (dim == 3) {
165         const CeedInt opt_elems       = block_size / (Q_1d * Q_1d);
166         const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1;
167         const CeedInt grid_size       = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0);
168 
169         CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args));
170       }
171     } break;
172     case CEED_EVAL_NONE: /* handled separately below */
173       break;
174     // LCOV_EXCL_START
175     case CEED_EVAL_DIV:
176     case CEED_EVAL_CURL:
177       return CeedError(ceed, CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]);
178       // LCOV_EXCL_STOP
179   }
180 
181   // Restore vectors, cover CEED_EVAL_NONE
182   CeedCallBackend(CeedVectorRestoreArray(v, &d_v));
183   if (eval_mode == CEED_EVAL_NONE) CeedCallBackend(CeedVectorSetArray(v, CEED_MEM_DEVICE, CEED_COPY_VALUES, (CeedScalar *)d_u));
184   if (eval_mode != CEED_EVAL_WEIGHT) CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u));
185   return CEED_ERROR_SUCCESS;
186 }
187 
188 //------------------------------------------------------------------------------
189 // Destroy basis
190 //------------------------------------------------------------------------------
191 static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) {
192   Ceed                   ceed;
193   CeedBasis_Cuda_shared *data;
194 
195   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
196   CeedCallBackend(CeedBasisGetData(basis, &data));
197   CeedCallCuda(ceed, cuModuleUnload(data->module));
198   CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d));
199   CeedCallCuda(ceed, cudaFree(data->d_interp_1d));
200   CeedCallCuda(ceed, cudaFree(data->d_grad_1d));
201   CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d));
202   CeedCallBackend(CeedFree(&data));
203   return CEED_ERROR_SUCCESS;
204 }
205 
206 //------------------------------------------------------------------------------
207 // Create tensor basis
208 //------------------------------------------------------------------------------
209 int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d,
210                                         const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) {
211   Ceed                   ceed;
212   char                  *basis_kernel_path, *basis_kernel_source;
213   CeedInt                num_comp;
214   const CeedInt          q_bytes      = Q_1d * sizeof(CeedScalar);
215   const CeedInt          interp_bytes = q_bytes * P_1d;
216   CeedBasis_Cuda_shared *data;
217 
218   CeedCallBackend(CeedBasisGetCeed(basis, &ceed));
219   CeedCallBackend(CeedCalloc(1, &data));
220 
221   // Copy basis data to GPU
222   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes));
223   CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice));
224   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes));
225   CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice));
226   CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes));
227   CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice));
228 
229   // Compute collocated gradient and copy to GPU
230   data->d_collo_grad_1d    = NULL;
231   bool has_collocated_grad = dim == 3 && Q_1d >= P_1d;
232 
233   if (has_collocated_grad) {
234     CeedScalar *collo_grad_1d;
235 
236     CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d));
237     CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d));
238     CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d));
239     CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice));
240     CeedCallBackend(CeedFree(&collo_grad_1d));
241   }
242 
243   // Compile basis kernels
244   CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp));
245   CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path));
246   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n");
247   CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source));
248   CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n");
249   CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D",
250                                    CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim),
251                                    "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad));
252   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp));
253   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose));
254   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad));
255   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose));
256   CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight));
257   CeedCallBackend(CeedFree(&basis_kernel_path));
258   CeedCallBackend(CeedFree(&basis_kernel_source));
259 
260   CeedCallBackend(CeedBasisSetData(basis, data));
261 
262   // Register backend functions
263   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared));
264   CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared));
265   return CEED_ERROR_SUCCESS;
266 }
267 
268 //------------------------------------------------------------------------------
269