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 CeedCudaInitInterp(CeedScalar *d_B, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B); 24 int CeedCudaInitGrad(CeedScalar *d_B, CeedScalar *d_G, CeedInt P_1d, CeedInt Q_1d, CeedScalar **c_B_ptr, CeedScalar **c_G_ptr); 25 int CeedCudaInitCollocatedGrad(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 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 34 Ceed_Cuda *ceed_Cuda; 35 CeedCallBackend(CeedGetData(ceed, &ceed_Cuda)); 36 CeedBasis_Cuda_shared *data; 37 CeedCallBackend(CeedBasisGetData(basis, &data)); 38 CeedInt dim, num_comp; 39 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 40 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 41 42 // Read vectors 43 const CeedScalar *d_u; 44 CeedScalar *d_v; 45 if (eval_mode != CEED_EVAL_WEIGHT) { 46 CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 47 } 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 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 55 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 56 CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 57 CeedCallBackend(CeedCudaInitInterp(data->d_interp_1d, P_1d, Q_1d, &data->c_B)); 58 void *interp_args[] = {(void *)&num_elem, &data->c_B, &d_u, &d_v}; 59 if (dim == 1) { 60 CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 61 1)); // avoid >512 total threads 62 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 63 CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 64 if (t_mode == CEED_TRANSPOSE) { 65 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 66 } else { 67 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 68 } 69 } else if (dim == 2) { 70 const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 71 // elems_per_block must be at least 1 72 CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 73 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 74 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 75 if (t_mode == CEED_TRANSPOSE) { 76 CeedCallBackend( 77 CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 78 } else { 79 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 80 } 81 } else if (dim == 3) { 82 CeedInt elems_per_block = 1; 83 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 84 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 85 if (t_mode == CEED_TRANSPOSE) { 86 CeedCallBackend( 87 CeedRunKernelDimSharedCuda(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 88 } else { 89 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 90 } 91 } 92 } break; 93 case CEED_EVAL_GRAD: { 94 CeedInt P_1d, Q_1d; 95 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 96 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 97 CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 98 if (data->d_collo_grad_1d) { 99 CeedCallBackend(CeedCudaInitCollocatedGrad(data->d_interp_1d, data->d_collo_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 100 } else { 101 CeedCallBackend(CeedCudaInitGrad(data->d_interp_1d, data->d_grad_1d, P_1d, Q_1d, &data->c_B, &data->c_G)); 102 } 103 void *grad_args[] = {(void *)&num_elem, &data->c_B, &data->c_G, &d_u, &d_v}; 104 if (dim == 1) { 105 CeedInt elems_per_block = CeedIntMin(ceed_Cuda->device_prop.maxThreadsDim[2], CeedIntMax(512 / thread_1d, 106 1)); // avoid >512 total threads 107 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 108 CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 109 if (t_mode == CEED_TRANSPOSE) { 110 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 111 } else { 112 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 113 } 114 } else if (dim == 2) { 115 const CeedInt opt_elems[7] = {0, 32, 8, 6, 4, 2, 8}; 116 // elems_per_block must be at least 1 117 CeedInt elems_per_block = CeedIntMax(thread_1d < 7 ? opt_elems[thread_1d] / num_comp : 1, 1); 118 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 119 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 120 if (t_mode == CEED_TRANSPOSE) { 121 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 122 } else { 123 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 124 } 125 } else if (dim == 3) { 126 CeedInt elems_per_block = 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 if (t_mode == CEED_TRANSPOSE) { 130 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 131 } else { 132 CeedCallBackend(CeedRunKernelDimSharedCuda(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 133 } 134 } 135 } break; 136 case CEED_EVAL_WEIGHT: { 137 CeedInt Q_1d; 138 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 139 void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 140 if (dim == 1) { 141 const CeedInt elems_per_block = 32 / Q_1d; 142 const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 143 CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, elems_per_block, 1, weight_args)); 144 } else if (dim == 2) { 145 const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 146 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 147 const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 148 CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args)); 149 } else if (dim == 3) { 150 const CeedInt opt_elems = 32 / (Q_1d * Q_1d); 151 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 152 const CeedInt gridsize = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 153 CeedCallBackend(CeedRunKernelDimCuda(ceed, data->Weight, gridsize, Q_1d, Q_1d, elems_per_block, weight_args)); 154 } 155 } break; 156 // LCOV_EXCL_START 157 // Evaluate the divergence to/from the quadrature points 158 case CEED_EVAL_DIV: 159 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported"); 160 // Evaluate the curl to/from the quadrature points 161 case CEED_EVAL_CURL: 162 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported"); 163 // Take no action, BasisApply should not have been called 164 case CEED_EVAL_NONE: 165 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 166 // LCOV_EXCL_STOP 167 } 168 169 // Restore vectors 170 if (eval_mode != CEED_EVAL_WEIGHT) { 171 CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 172 } 173 CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 174 return CEED_ERROR_SUCCESS; 175 } 176 177 //------------------------------------------------------------------------------ 178 // Destroy basis 179 //------------------------------------------------------------------------------ 180 static int CeedBasisDestroy_Cuda_shared(CeedBasis basis) { 181 Ceed ceed; 182 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 183 184 CeedBasis_Cuda_shared *data; 185 CeedCallBackend(CeedBasisGetData(basis, &data)); 186 187 CeedCallCuda(ceed, cuModuleUnload(data->module)); 188 189 CeedCallCuda(ceed, cudaFree(data->d_q_weight_1d)); 190 CeedCallCuda(ceed, cudaFree(data->d_interp_1d)); 191 CeedCallCuda(ceed, cudaFree(data->d_grad_1d)); 192 CeedCallCuda(ceed, cudaFree(data->d_collo_grad_1d)); 193 194 CeedCallBackend(CeedFree(&data)); 195 196 return CEED_ERROR_SUCCESS; 197 } 198 199 //------------------------------------------------------------------------------ 200 // Create tensor basis 201 //------------------------------------------------------------------------------ 202 int CeedBasisCreateTensorH1_Cuda_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 203 const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 204 Ceed ceed; 205 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 206 CeedBasis_Cuda_shared *data; 207 CeedCallBackend(CeedCalloc(1, &data)); 208 209 // Copy basis data to GPU 210 const CeedInt q_bytes = Q_1d * sizeof(CeedScalar); 211 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_q_weight_1d, q_bytes)); 212 CeedCallCuda(ceed, cudaMemcpy(data->d_q_weight_1d, q_weight_1d, q_bytes, cudaMemcpyHostToDevice)); 213 214 const CeedInt interp_bytes = q_bytes * P_1d; 215 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_interp_1d, interp_bytes)); 216 CeedCallCuda(ceed, cudaMemcpy(data->d_interp_1d, interp_1d, interp_bytes, cudaMemcpyHostToDevice)); 217 218 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_grad_1d, interp_bytes)); 219 CeedCallCuda(ceed, cudaMemcpy(data->d_grad_1d, grad_1d, interp_bytes, cudaMemcpyHostToDevice)); 220 221 // Compute collocated gradient and copy to GPU 222 data->d_collo_grad_1d = NULL; 223 bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 224 if (has_collocated_grad) { 225 CeedScalar *collo_grad_1d; 226 CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 227 CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 228 CeedCallCuda(ceed, cudaMalloc((void **)&data->d_collo_grad_1d, q_bytes * Q_1d)); 229 CeedCallCuda(ceed, cudaMemcpy(data->d_collo_grad_1d, collo_grad_1d, q_bytes * Q_1d, cudaMemcpyHostToDevice)); 230 CeedCallBackend(CeedFree(&collo_grad_1d)); 231 } 232 233 // Compile basis kernels 234 CeedInt num_comp; 235 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 236 char *basis_kernel_path, *basis_kernel_source; 237 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor.h", &basis_kernel_path)); 238 CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source -----\n"); 239 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 240 CeedDebug256(ceed, 2, "----- Loading Basis Kernel Source Complete -----\n"); 241 CeedCallBackend(CeedCompileCuda(ceed, basis_kernel_source, &data->module, 8, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", CeedIntMax(Q_1d, P_1d), 242 "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), "BASIS_NUM_QPTS", 243 CeedIntPow(Q_1d, dim), "BASIS_HAS_COLLOCATED_GRAD", has_collocated_grad)); 244 CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Interp", &data->Interp)); 245 CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 246 CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Grad", &data->Grad)); 247 CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "GradTranspose", &data->GradTranspose)); 248 CeedCallBackend(CeedGetKernelCuda(ceed, data->module, "Weight", &data->Weight)); 249 CeedCallBackend(CeedFree(&basis_kernel_path)); 250 CeedCallBackend(CeedFree(&basis_kernel_source)); 251 252 CeedCallBackend(CeedBasisSetData(basis, data)); 253 254 // Register backend functions 255 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Cuda_shared)); 256 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Cuda_shared)); 257 return CEED_ERROR_SUCCESS; 258 } 259 260 //------------------------------------------------------------------------------ 261