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 // Read vectors 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 150 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 151 void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 152 if (dim == 1) { 153 const CeedInt elems_per_block = 32 / Q_1d; 154 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 155 156 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args)); 157 } else if (dim == 2) { 158 const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 159 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 160 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 161 162 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 163 } else if (dim == 3) { 164 const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 165 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 166 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 167 168 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 169 } 170 } break; 171 // LCOV_EXCL_START 172 // Evaluate the divergence to/from the quadrature points 173 case CEED_EVAL_DIV: 174 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported"); 175 // Evaluate the curl to/from the quadrature points 176 case CEED_EVAL_CURL: 177 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported"); 178 // Take no action, BasisApply should not have been called 179 case CEED_EVAL_NONE: 180 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 181 // LCOV_EXCL_STOP 182 } 183 184 // Restore vectors 185 if (eval_mode != CEED_EVAL_WEIGHT) { 186 CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 187 } 188 CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 189 return CEED_ERROR_SUCCESS; 190 } 191 192 //------------------------------------------------------------------------------ 193 // Destroy basis 194 //------------------------------------------------------------------------------ 195 static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 196 Ceed ceed; 197 CeedBasis_Cuda_shared *data; 198 199 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 200 CeedCallBackend(CeedBasisGetData(basis, &data)); 201 CeedCallCuda(ceed, cuModuleUnload(data->module)); 202 CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 203 CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 204 CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 205 CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 206 CeedCallBackend(CeedFree(&data)); 207 return CEED_ERROR_SUCCESS; 208 } 209 210 //------------------------------------------------------------------------------ 211 // Create tensor basis 212 //------------------------------------------------------------------------------ 213 int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 214 const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 215 Ceed ceed; 216 char *basis_kernel_path, *basis_kernel_source; 217 CeedInt num_comp; 218 const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 219 const CeedInt interp_bytes = q_bytes * P_1d; 220 CeedBasis_Cuda_shared *data; 221 222 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 223 CeedCallBackend(CeedCalloc(1, &data)); 224 225 // Copy basis data to GPU 226 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 227 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 228 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 229 CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 230 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 231 CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 232 233 // Compute collocated gradient and copy to GPU 234 data->d_collo_grad_1d = NULL; 235 bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 236 237 if (has_collocated_grad) { 238 CeedScalar *collo_grad_1d; 239 240 CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 241 CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 242 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 243 CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 244 CeedCallBackend(CeedFree(&collo_grad_1d)); 245 } 246 247 // Compile basis kernels 248 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 249 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path)); 250 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 251 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 252 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete -----\n"); 253 CeedCallBackend(CeedCompile_Cuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 254 CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 255 "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 256 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Interp", &data->Interp)); 257 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 258 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Grad", &data->Grad)); 259 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 260 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, "Weight", &data->Weight)); 261 CeedCallBackend(CeedFree(&basis_kernel_path)); 262 CeedCallBackend(CeedFree(&basis_kernel_source)); 263 264 CeedCallBackend(CeedBasisSetData(basis, data)); 265 266 // Register backend functions 267 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 268 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 269 return CEED_ERROR_SUCCESS; 270 } 271 272 //------------------------------------------------------------------------------ 273