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 <stdbool.h> 12 #include <stddef.h> 13 #include <hip/hip_runtime.h> 14 15 #include "../hip/ceed-hip-common.h" 16 #include "../hip/ceed-hip-compile.h" 17 #include "ceed-hip-shared.h" 18 19 //------------------------------------------------------------------------------ 20 // Compute a block size based on required minimum threads 21 //------------------------------------------------------------------------------ 22 static CeedInt ComputeBlockSizeFromRequirement(const CeedInt required) { 23 CeedInt maxSize = 1024; // Max total threads per block 24 CeedInt currentSize = 64; // Start with one group 25 26 while (currentSize < maxSize) { 27 if (currentSize > required) break; 28 else currentSize = currentSize * 2; 29 } 30 return currentSize; 31 } 32 33 //------------------------------------------------------------------------------ 34 // Compute required thread block sizes for basis kernels given P, Q, dim, and 35 // num_comp (num_comp not currently used, but may be again in other basis 36 // parallelization options) 37 //------------------------------------------------------------------------------ 38 static int ComputeBasisThreadBlockSizes(const CeedInt dim, const CeedInt P_1d, const CeedInt Q_1d, const CeedInt num_comp, CeedInt *block_sizes) { 39 // Note that this will use the same block sizes for all dimensions when compiling, 40 // but as each basis object is defined for a particular dimension, we will never 41 // call any kernels except the ones for the dimension for which we have computed the 42 // block sizes. 43 const CeedInt thread_1d = CeedIntMax(P_1d, Q_1d); 44 switch (dim) { 45 case 1: { 46 // Interp kernels: 47 block_sizes[0] = 256; 48 49 // Grad kernels: 50 block_sizes[1] = 256; 51 52 // Weight kernels: 53 block_sizes[2] = 256; 54 } break; 55 case 2: { 56 // Interp kernels: 57 CeedInt required = thread_1d * thread_1d; 58 block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 59 60 // Grad kernels: currently use same required minimum threads 61 block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 62 63 // Weight kernels: 64 required = CeedIntMax(64, Q_1d * Q_1d); 65 block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 66 67 } break; 68 case 3: { 69 // Interp kernels: 70 CeedInt required = thread_1d * thread_1d; 71 block_sizes[0] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 72 73 // Grad kernels: currently use same required minimum threads 74 block_sizes[1] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 75 76 // Weight kernels: 77 required = Q_1d * Q_1d * Q_1d; 78 block_sizes[2] = CeedIntMax(256, ComputeBlockSizeFromRequirement(required)); 79 } 80 } 81 82 return CEED_ERROR_SUCCESS; 83 } 84 85 //------------------------------------------------------------------------------ 86 // Apply basis 87 //------------------------------------------------------------------------------ 88 int CeedBasisApplyTensor_Hip_shared(CeedBasis basis, const CeedInt num_elem, CeedTransposeMode t_mode, CeedEvalMode eval_mode, CeedVector u, 89 CeedVector v) { 90 Ceed ceed; 91 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 92 Ceed_Hip *ceed_Hip; 93 CeedCallBackend(CeedGetData(ceed, &ceed_Hip)); 94 CeedBasis_Hip_shared *data; 95 CeedCallBackend(CeedBasisGetData(basis, &data)); 96 CeedInt dim, num_comp; 97 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 98 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 99 100 // Read vectors 101 const CeedScalar *d_u; 102 CeedScalar *d_v; 103 if (u != CEED_VECTOR_NONE) CeedCallBackend(CeedVectorGetArrayRead(u, CEED_MEM_DEVICE, &d_u)); 104 else CeedCheck(eval_mode == CEED_EVAL_WEIGHT, ceed, CEED_ERROR_BACKEND, "An input vector is required for this CeedEvalMode"); 105 CeedCallBackend(CeedVectorGetArrayWrite(v, CEED_MEM_DEVICE, &d_v)); 106 107 // Apply basis operation 108 switch (eval_mode) { 109 case CEED_EVAL_INTERP: { 110 CeedInt P_1d, Q_1d; 111 CeedInt block_size = data->block_sizes[0]; 112 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 113 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 114 CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 115 void *interp_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_u, &d_v}; 116 if (dim == 1) { 117 CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64; 118 elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 119 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 120 CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 121 122 if (t_mode == CEED_TRANSPOSE) { 123 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 124 } else { 125 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, 1, elems_per_block, shared_mem, interp_args)); 126 } 127 } else if (dim == 2) { 128 // Check if required threads is small enough to do multiple elems 129 const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1); 130 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 131 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 132 133 if (t_mode == CEED_TRANSPOSE) { 134 CeedCallBackend( 135 CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 136 } else { 137 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 138 } 139 } else if (dim == 3) { 140 const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1); 141 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 142 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 143 144 if (t_mode == CEED_TRANSPOSE) { 145 CeedCallBackend( 146 CeedRunKernelDimShared_Hip(ceed, data->InterpTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 147 } else { 148 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Interp, grid, thread_1d, thread_1d, elems_per_block, shared_mem, interp_args)); 149 } 150 } 151 } break; 152 case CEED_EVAL_GRAD: { 153 CeedInt P_1d, Q_1d; 154 CeedInt block_size = data->block_sizes[1]; 155 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 156 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 157 CeedInt thread_1d = CeedIntMax(Q_1d, P_1d); 158 CeedScalar *d_grad_1d = data->d_grad_1d; 159 if (data->d_collo_grad_1d) { 160 d_grad_1d = data->d_collo_grad_1d; 161 } 162 void *grad_args[] = {(void *)&num_elem, &data->d_interp_1d, &d_grad_1d, &d_u, &d_v}; 163 if (dim == 1) { 164 CeedInt elems_per_block = 64 * thread_1d > 256 ? 256 / thread_1d : 64; 165 elems_per_block = elems_per_block > 0 ? elems_per_block : 1; 166 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 167 CeedInt shared_mem = elems_per_block * thread_1d * sizeof(CeedScalar); 168 169 if (t_mode == CEED_TRANSPOSE) { 170 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 171 } else { 172 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, 1, elems_per_block, shared_mem, grad_args)); 173 } 174 } else if (dim == 2) { 175 // Check if required threads is small enough to do multiple elems 176 const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1); 177 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 178 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 179 180 if (t_mode == CEED_TRANSPOSE) { 181 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 182 } else { 183 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 184 } 185 } else if (dim == 3) { 186 const CeedInt elems_per_block = CeedIntMax(block_size / (thread_1d * thread_1d), 1); 187 CeedInt grid = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 188 CeedInt shared_mem = elems_per_block * thread_1d * thread_1d * sizeof(CeedScalar); 189 190 if (t_mode == CEED_TRANSPOSE) { 191 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->GradTranspose, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 192 } else { 193 CeedCallBackend(CeedRunKernelDimShared_Hip(ceed, data->Grad, grid, thread_1d, thread_1d, elems_per_block, shared_mem, grad_args)); 194 } 195 } 196 } break; 197 case CEED_EVAL_WEIGHT: { 198 CeedInt Q_1d; 199 CeedInt block_size = data->block_sizes[2]; 200 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 201 void *weight_args[] = {(void *)&num_elem, (void *)&data->d_q_weight_1d, &d_v}; 202 if (dim == 1) { 203 const CeedInt opt_elems = block_size / Q_1d; 204 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 205 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 206 207 CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, elems_per_block, 1, weight_args)); 208 } else if (dim == 2) { 209 const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 210 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 211 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 212 213 CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 214 } else if (dim == 3) { 215 const CeedInt opt_elems = block_size / (Q_1d * Q_1d); 216 const CeedInt elems_per_block = opt_elems > 0 ? opt_elems : 1; 217 const CeedInt grid_size = num_elem / elems_per_block + ((num_elem / elems_per_block * elems_per_block < num_elem) ? 1 : 0); 218 219 CeedCallBackend(CeedRunKernelDim_Hip(ceed, data->Weight, grid_size, Q_1d, Q_1d, elems_per_block, weight_args)); 220 } 221 } break; 222 // LCOV_EXCL_START 223 // Evaluate the divergence to/from the quadrature points 224 case CEED_EVAL_DIV: 225 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_DIV not supported"); 226 // Evaluate the curl to/from the quadrature points 227 case CEED_EVAL_CURL: 228 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_CURL not supported"); 229 // Take no action, BasisApply should not have been called 230 case CEED_EVAL_NONE: 231 return CeedError(ceed, CEED_ERROR_BACKEND, "CEED_EVAL_NONE does not make sense in this context"); 232 // LCOV_EXCL_STOP 233 } 234 235 // Restore vectors 236 if (eval_mode != CEED_EVAL_WEIGHT) { 237 CeedCallBackend(CeedVectorRestoreArrayRead(u, &d_u)); 238 } 239 CeedCallBackend(CeedVectorRestoreArray(v, &d_v)); 240 return CEED_ERROR_SUCCESS; 241 } 242 243 //------------------------------------------------------------------------------ 244 // Destroy basis 245 //------------------------------------------------------------------------------ 246 static int CeedBasisDestroy_Hip_shared(CeedBasis basis) { 247 Ceed ceed; 248 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 249 250 CeedBasis_Hip_shared *data; 251 CeedCallBackend(CeedBasisGetData(basis, &data)); 252 253 CeedCallHip(ceed, hipModuleUnload(data->module)); 254 255 CeedCallHip(ceed, hipFree(data->d_q_weight_1d)); 256 CeedCallHip(ceed, hipFree(data->d_interp_1d)); 257 CeedCallHip(ceed, hipFree(data->d_grad_1d)); 258 CeedCallHip(ceed, hipFree(data->d_collo_grad_1d)); 259 CeedCallBackend(CeedFree(&data)); 260 261 return CEED_ERROR_SUCCESS; 262 } 263 264 //------------------------------------------------------------------------------ 265 // Create tensor basis 266 //------------------------------------------------------------------------------ 267 int CeedBasisCreateTensorH1_Hip_shared(CeedInt dim, CeedInt P_1d, CeedInt Q_1d, const CeedScalar *interp_1d, const CeedScalar *grad_1d, 268 const CeedScalar *q_ref_1d, const CeedScalar *q_weight_1d, CeedBasis basis) { 269 Ceed ceed; 270 CeedCallBackend(CeedBasisGetCeed(basis, &ceed)); 271 CeedBasis_Hip_shared *data; 272 CeedCallBackend(CeedCalloc(1, &data)); 273 274 // Copy basis data to GPU 275 const CeedInt qBytes = Q_1d * sizeof(CeedScalar); 276 CeedCallHip(ceed, hipMalloc((void **)&data->d_q_weight_1d, qBytes)); 277 CeedCallHip(ceed, hipMemcpy(data->d_q_weight_1d, q_weight_1d, qBytes, hipMemcpyHostToDevice)); 278 279 const CeedInt iBytes = qBytes * P_1d; 280 CeedCallHip(ceed, hipMalloc((void **)&data->d_interp_1d, iBytes)); 281 CeedCallHip(ceed, hipMemcpy(data->d_interp_1d, interp_1d, iBytes, hipMemcpyHostToDevice)); 282 283 CeedCallHip(ceed, hipMalloc((void **)&data->d_grad_1d, iBytes)); 284 CeedCallHip(ceed, hipMemcpy(data->d_grad_1d, grad_1d, iBytes, hipMemcpyHostToDevice)); 285 286 // Compute collocated gradient and copy to GPU 287 data->d_collo_grad_1d = NULL; 288 bool has_collocated_grad = dim == 3 && Q_1d >= P_1d; 289 if (has_collocated_grad) { 290 CeedScalar *collo_grad_1d; 291 CeedCallBackend(CeedMalloc(Q_1d * Q_1d, &collo_grad_1d)); 292 CeedCallBackend(CeedBasisGetCollocatedGrad(basis, collo_grad_1d)); 293 CeedCallHip(ceed, hipMalloc((void **)&data->d_collo_grad_1d, qBytes * Q_1d)); 294 CeedCallHip(ceed, hipMemcpy(data->d_collo_grad_1d, collo_grad_1d, qBytes * Q_1d, hipMemcpyHostToDevice)); 295 CeedCallBackend(CeedFree(&collo_grad_1d)); 296 } 297 298 // Set number of threads per block for basis kernels 299 CeedInt num_comp; 300 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 301 CeedCallBackend(ComputeBasisThreadBlockSizes(dim, P_1d, Q_1d, num_comp, data->block_sizes)); 302 303 // Compile basis kernels 304 char *basis_kernel_path, *basis_kernel_source; 305 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/hip/hip-shared-basis-tensor.h", &basis_kernel_path)); 306 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source -----\n"); 307 CeedCallBackend(CeedLoadSourceToBuffer(ceed, basis_kernel_path, &basis_kernel_source)); 308 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Basis Kernel Source Complete! -----\n"); 309 CeedCallBackend(CeedCompile_Hip(ceed, basis_kernel_source, &data->module, 11, "BASIS_Q_1D", Q_1d, "BASIS_P_1D", P_1d, "T_1D", 310 CeedIntMax(Q_1d, P_1d), "BASIS_DIM", dim, "BASIS_NUM_COMP", num_comp, "BASIS_NUM_NODES", CeedIntPow(P_1d, dim), 311 "BASIS_NUM_QPTS", CeedIntPow(Q_1d, dim), "BASIS_INTERP_BLOCK_SIZE", data->block_sizes[0], "BASIS_GRAD_BLOCK_SIZE", 312 data->block_sizes[1], "BASIS_WEIGHT_BLOCK_SIZE", data->block_sizes[2], "BASIS_HAS_COLLOCATED_GRAD", 313 has_collocated_grad)); 314 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Interp", &data->Interp)); 315 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "InterpTranspose", &data->InterpTranspose)); 316 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Grad", &data->Grad)); 317 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "GradTranspose", &data->GradTranspose)); 318 CeedCallBackend(CeedGetKernel_Hip(ceed, data->module, "Weight", &data->Weight)); 319 CeedCallBackend(CeedFree(&basis_kernel_path)); 320 CeedCallBackend(CeedFree(&basis_kernel_source)); 321 322 CeedCallBackend(CeedBasisSetData(basis, data)); 323 324 // Register backend functions 325 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Apply", CeedBasisApplyTensor_Hip_shared)); 326 CeedCallBackend(CeedSetBackendFunction(ceed, "Basis", basis, "Destroy", CeedBasisDestroy_Hip_shared)); 327 return CEED_ERROR_SUCCESS; 328 } 329 330 //------------------------------------------------------------------------------ 331