1 // Copyright (c) 2017-2024, 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 #define CEED_DEBUG_COLOR 12 9 10 #include <ceed.h> 11 #include <ceed/backend.h> 12 #include <ceed/jit-tools.h> 13 #include <cuda_runtime.h> 14 15 #include <iostream> 16 #include <sstream> 17 #include <string> 18 19 #include "../cuda-ref/ceed-cuda-ref.h" 20 #include "../cuda-shared/ceed-cuda-shared.h" 21 #include "../cuda/ceed-cuda-common.h" 22 #include "../cuda/ceed-cuda-compile.h" 23 #include "ceed-cuda-gen.h" 24 25 //------------------------------------------------------------------------------ 26 // Build single operator kernel 27 //------------------------------------------------------------------------------ 28 extern "C" int CeedOperatorBuildKernel_Cuda_gen(CeedOperator op) { 29 using std::ostringstream; 30 using std::string; 31 32 bool is_setup_done, is_identity_qf; 33 struct cudaDeviceProp prop; 34 Ceed ceed; 35 Ceed_Cuda *ceed_data; 36 CeedSize l_size; 37 CeedInt Q, P_1d = 0, Q_1d = 0, elem_size, num_input_fields, num_output_fields, num_comp, dim = 1; 38 CeedEvalMode eval_mode; 39 CeedElemRestriction elem_rstr; 40 CeedElemRestriction_Cuda *rstr_data; 41 CeedBasis basis; 42 CeedBasis_Cuda_shared *basis_data; 43 CeedQFunctionField *qf_input_fields, *qf_output_fields; 44 CeedQFunction_Cuda_gen *qf_data; 45 CeedQFunction qf; 46 CeedOperatorField *op_input_fields, *op_output_fields; 47 CeedOperator_Cuda_gen *data; 48 49 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 50 if (is_setup_done) return CEED_ERROR_SUCCESS; 51 52 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 53 CeedCallBackend(CeedOperatorGetData(op, &data)); 54 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 55 CeedCallBackend(CeedQFunctionGetData(qf, &qf_data)); 56 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 57 Q_1d = Q; 58 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 59 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 60 61 // TODO: put in a function? 62 // Check for restriction only identity operator 63 CeedCallBackend(CeedQFunctionIsIdentity(qf, &is_identity_qf)); 64 if (is_identity_qf) { 65 CeedEvalMode eval_mode_in, eval_mode_out; 66 67 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[0], &eval_mode_in)); 68 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[0], &eval_mode_out)); 69 CeedCheck(eval_mode_in != CEED_EVAL_NONE || eval_mode_out != CEED_EVAL_NONE, ceed, CEED_ERROR_BACKEND, 70 "Backend does not implement restriction only identity operators"); 71 } 72 73 ostringstream code; 74 75 // TODO: put in a function? 76 // Add atomicAdd function for old NVidia architectures 77 CeedCallBackend(CeedGetData(ceed, &ceed_data)); 78 CeedCallBackend(cudaGetDeviceProperties(&prop, ceed_data->device_id)); 79 if ((prop.major < 6) && (CEED_SCALAR_TYPE != CEED_SCALAR_FP32)) { 80 char *atomic_add_source; 81 const char *atomic_add_path; 82 83 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-atomic-add-fallback.h", &atomic_add_path)); 84 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Atomic Add Source -----\n"); 85 CeedCallBackend(CeedLoadSourceToBuffer(ceed, atomic_add_path, &atomic_add_source)); 86 code << atomic_add_source; 87 CeedCallBackend(CeedFree(&atomic_add_path)); 88 CeedCallBackend(CeedFree(&atomic_add_source)); 89 } 90 91 // Load basis source files 92 // TODO: generalize to accept different device functions? 93 { 94 char *tensor_basis_kernel_source; 95 const char *tensor_basis_kernel_path; 96 97 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-shared-basis-tensor-templates.h", &tensor_basis_kernel_path)); 98 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Tensor Basis Kernel Source -----\n"); 99 CeedCallBackend(CeedLoadSourceToBuffer(ceed, tensor_basis_kernel_path, &tensor_basis_kernel_source)); 100 code << tensor_basis_kernel_source; 101 CeedCallBackend(CeedFree(&tensor_basis_kernel_path)); 102 CeedCallBackend(CeedFree(&tensor_basis_kernel_source)); 103 } 104 { 105 char *cuda_gen_template_source; 106 const char *cuda_gen_template_path; 107 108 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-gen-templates.h", &cuda_gen_template_path)); 109 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Cuda-Gen Template Source -----\n"); 110 CeedCallBackend(CeedLoadSourceToBuffer(ceed, cuda_gen_template_path, &cuda_gen_template_source)); 111 code << cuda_gen_template_source; 112 CeedCallBackend(CeedFree(&cuda_gen_template_path)); 113 CeedCallBackend(CeedFree(&cuda_gen_template_source)); 114 } 115 116 // Get QFunction source and name 117 string qfunction_source(qf_data->qfunction_source); 118 string qfunction_name(qf_data->qfunction_name); 119 string operator_name; 120 operator_name = "CeedKernelCudaGenOperator_" + qfunction_name; 121 122 // Find dim, P_1d, Q_1d 123 data->max_P_1d = 0; 124 for (CeedInt i = 0; i < num_input_fields; i++) { 125 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 126 if (basis != CEED_BASIS_NONE) { 127 bool is_tensor; 128 129 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 130 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 131 132 // Collect dim, P_1d, and Q_1d 133 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 134 CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 135 CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 136 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 137 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 138 data->max_P_1d = CeedIntMax(data->max_P_1d, P_1d); 139 } 140 } 141 // Check output bases for Q_1d, dim as well 142 // The only input basis might be CEED_BASIS_NONE 143 for (CeedInt i = 0; i < num_output_fields; i++) { 144 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 145 if (basis != CEED_BASIS_NONE) { 146 bool is_tensor; 147 148 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 149 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 150 151 // Collect Q_1d 152 CeedCallBackend(CeedBasisGetDimension(basis, &dim)); 153 CeedCallBackend(CeedBasisIsTensor(basis, &is_tensor)); 154 CeedCheck(is_tensor, ceed, CEED_ERROR_BACKEND, "Backend does not implement operators with non-tensor basis"); 155 CeedCallBackend(CeedBasisGetNumQuadraturePoints1D(basis, &Q_1d)); 156 } 157 } 158 data->dim = dim; 159 data->Q_1d = Q_1d; 160 161 // Only use 3D collocated gradient parallelization strategy when gradient is computed 162 // TODO: put in a function? 163 bool use_collograd_parallelization = false; 164 165 if (dim == 3) { 166 bool was_grad_found = false; 167 168 for (CeedInt i = 0; i < num_input_fields; i++) { 169 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 170 if (eval_mode == CEED_EVAL_GRAD) { 171 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 172 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 173 use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true); 174 was_grad_found = true; 175 } 176 } 177 for (CeedInt i = 0; i < num_output_fields; i++) { 178 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 179 if (eval_mode == CEED_EVAL_GRAD) { 180 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 181 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 182 use_collograd_parallelization = basis_data->d_collo_grad_1d && (was_grad_found ? use_collograd_parallelization : true); 183 was_grad_found = true; 184 } 185 } 186 } 187 188 // Define CEED_Q_VLA 189 code << "\n#undef CEED_Q_VLA\n"; 190 if (dim != 3 || use_collograd_parallelization) { 191 code << "#define CEED_Q_VLA 1\n\n"; 192 } else { 193 code << "#define CEED_Q_VLA " << Q_1d << "\n\n"; 194 } 195 196 code << qfunction_source; 197 198 // Setup 199 code << "\n// -----------------------------------------------------------------------------\n"; 200 code << "\nextern \"C\" __global__ void " << operator_name 201 << "(CeedInt num_elem, void* ctx, FieldsInt_Cuda indices, Fields_Cuda fields, Fields_Cuda B, Fields_Cuda G, CeedScalar* W) {\n"; 202 for (CeedInt i = 0; i < num_input_fields; i++) { 203 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 204 if (eval_mode != CEED_EVAL_WEIGHT) { // Skip CEED_EVAL_WEIGHT 205 code << " const CeedScalar* d_u_" << i << " = fields.inputs[" << i << "];\n"; 206 } 207 } 208 209 for (CeedInt i = 0; i < num_output_fields; i++) { 210 code << " CeedScalar* d_v_" << i << " = fields.outputs[" << i << "];\n"; 211 } 212 213 code << " const CeedInt dim = " << dim << ";\n"; 214 code << " const CeedInt Q_1d = " << Q_1d << ";\n"; 215 216 code << " extern __shared__ CeedScalar slice[];\n"; 217 // TODO put in a function? InitSharedData_Cuda? 218 code << " SharedData_Cuda data;\n"; 219 code << " data.t_id_x = threadIdx.x;\n"; 220 code << " data.t_id_y = threadIdx.y;\n"; 221 code << " data.t_id_z = threadIdx.z;\n"; 222 code << " data.t_id = threadIdx.x + threadIdx.y*blockDim.x + threadIdx.z*blockDim.y*blockDim.x;\n"; 223 code << " data.slice = slice+data.t_id_z*T_1D" << (dim > 1 ? "*T_1D" : "") << ";\n"; 224 225 code << "\n // -- Input field constants and basis data --\n"; 226 // TODO: Put in a function? 227 // Initialize constants, and matrices B and G 228 for (CeedInt i = 0; i < num_input_fields; i++) { 229 code << " // ---- Input field " << i << " ----\n"; 230 // Get elem_size, eval_mode, num_comp 231 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 232 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 233 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 234 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 235 236 // Set field constants 237 if (eval_mode != CEED_EVAL_WEIGHT) { 238 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 239 if (basis != CEED_BASIS_NONE) { 240 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 241 code << " const CeedInt P_in_" << i << " = " << P_1d << ";\n"; 242 } else { 243 code << " const CeedInt P_in_" << i << " = " << Q_1d << ";\n"; 244 } 245 code << " const CeedInt num_comp_in_" << i << " = " << num_comp << ";\n"; 246 } 247 248 // Load basis data 249 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 250 switch (eval_mode) { 251 case CEED_EVAL_NONE: 252 break; 253 case CEED_EVAL_INTERP: 254 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 255 data->B.inputs[i] = basis_data->d_interp_1d; 256 code << " __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n"; 257 code << " loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n"; 258 break; 259 case CEED_EVAL_GRAD: 260 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 261 data->B.inputs[i] = basis_data->d_interp_1d; 262 code << " __shared__ CeedScalar s_B_in_" << i << "[" << P_1d * Q_1d << "];\n"; 263 code << " loadMatrix<P_in_" << i << ",Q_1d>(data, B.inputs[" << i << "], s_B_in_" << i << ");\n"; 264 if (use_collograd_parallelization) { 265 data->G.inputs[i] = basis_data->d_collo_grad_1d; 266 code << " __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * Q_1d << "];\n"; 267 code << " loadMatrix<Q_1d,Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i << ");\n"; 268 } else { 269 bool has_collo_grad = basis_data->d_collo_grad_1d; 270 data->G.inputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 271 code << " __shared__ CeedScalar s_G_in_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n"; 272 code << " loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_in_" + std::to_string(i))) << ",Q_1d>(data, G.inputs[" << i << "], s_G_in_" << i 273 << ");\n"; 274 } 275 break; 276 case CEED_EVAL_WEIGHT: 277 break; // No action 278 case CEED_EVAL_DIV: 279 break; // TODO: Not implemented 280 case CEED_EVAL_CURL: 281 break; // TODO: Not implemented 282 } 283 } 284 285 code << "\n // -- Output field constants and basis data --\n"; 286 for (CeedInt i = 0; i < num_output_fields; i++) { 287 code << " // ---- Output field " << i << " ----\n"; 288 // Get elem_size, eval_mode, num_comp 289 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 290 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 291 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 292 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 293 294 // Set field constants 295 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 296 if (basis != CEED_BASIS_NONE) { 297 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 298 code << " const CeedInt P_out_" << i << " = " << P_1d << ";\n"; 299 } else { 300 code << " const CeedInt P_out_" << i << " = " << Q_1d << ";\n"; 301 } 302 code << " const CeedInt num_comp_out_" << i << " = " << num_comp << ";\n"; 303 304 // Load basis data 305 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 306 switch (eval_mode) { 307 case CEED_EVAL_NONE: 308 break; // No action 309 case CEED_EVAL_INTERP: 310 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 311 data->B.outputs[i] = basis_data->d_interp_1d; 312 code << " __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n"; 313 code << " loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n"; 314 break; 315 case CEED_EVAL_GRAD: 316 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 317 data->B.outputs[i] = basis_data->d_interp_1d; 318 code << " __shared__ CeedScalar s_B_out_" << i << "[" << P_1d * Q_1d << "];\n"; 319 code << " loadMatrix<P_out_" << i << ",Q_1d>(data, B.outputs[" << i << "], s_B_out_" << i << ");\n"; 320 if (use_collograd_parallelization) { 321 data->G.outputs[i] = basis_data->d_collo_grad_1d; 322 code << " __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * Q_1d << "];\n"; 323 code << " loadMatrix<Q_1d,Q_1d>(data, G.outputs[" << i << "], s_G_out_" << i << ");\n"; 324 } else { 325 bool has_collo_grad = basis_data->d_collo_grad_1d; 326 data->G.outputs[i] = has_collo_grad ? basis_data->d_collo_grad_1d : basis_data->d_grad_1d; 327 code << " __shared__ CeedScalar s_G_out_" << i << "[" << Q_1d * (has_collo_grad ? Q_1d : P_1d) << "];\n"; 328 code << " loadMatrix<" << (has_collo_grad ? "Q_1d" : ("P_out_" + std::to_string(i))) << ",Q_1d>(data, G.outputs[" << i << "], s_G_out_" 329 << i << ");\n"; 330 } 331 break; 332 // LCOV_EXCL_START 333 case CEED_EVAL_WEIGHT: { 334 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 335 break; // Should not occur 336 } 337 case CEED_EVAL_DIV: 338 case CEED_EVAL_CURL: { 339 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 340 break; // Should not occur 341 } 342 // LCOV_EXCL_STOP 343 } 344 } 345 code << "\n // -- Element loop --\n"; 346 code << " __syncthreads();\n"; 347 code << " for (CeedInt elem = blockIdx.x*blockDim.z + threadIdx.z; elem < num_elem; elem += gridDim.x*blockDim.z) {\n"; 348 // Input basis apply if needed 349 // Generate the correct eval mode code for each input 350 code << " // -- Input field restrictions and basis actions --\n"; 351 for (CeedInt i = 0; i < num_input_fields; i++) { 352 code << " // ---- Input field " << i << " ----\n"; 353 // Get elem_size, eval_mode, num_comp 354 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 355 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 356 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 357 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 358 359 // TODO: put in a function? 360 // Restriction 361 if (eval_mode != CEED_EVAL_WEIGHT && !((eval_mode == CEED_EVAL_NONE) && use_collograd_parallelization)) { 362 code << " CeedScalar r_u_" << i << "[num_comp_in_" << i << "*P_in_" << i << "];\n"; 363 364 bool is_strided; 365 366 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 367 if (!is_strided) { 368 CeedInt comp_stride; 369 370 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 371 code << " const CeedInt l_size_in_" << i << " = " << l_size << ";\n"; 372 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 373 code << " // CompStride: " << comp_stride << "\n"; 374 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 375 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 376 code << " readDofsOffset" << dim << "d<num_comp_in_" << i << ", " << comp_stride << ", P_in_" << i << ">(data, l_size_in_" << i 377 << ", elem, indices.inputs[" << i << "], d_u_" << i << ", r_u_" << i << ");\n"; 378 } else { 379 bool has_backend_strides; 380 CeedInt num_elem; 381 382 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 383 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 384 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 385 386 if (!has_backend_strides) { 387 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 388 } 389 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 390 code << " readDofsStrided" << dim << "d<num_comp_in_" << i << ",P_in_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2] 391 << ">(data, elem, d_u_" << i << ", r_u_" << i << ");\n"; 392 } 393 } 394 395 // TODO: put in a function? 396 // Basis action 397 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 398 switch (eval_mode) { 399 case CEED_EVAL_NONE: 400 if (!use_collograd_parallelization) { 401 code << " CeedScalar* r_t_" << i << " = r_u_" << i << ";\n"; 402 } 403 break; 404 case CEED_EVAL_INTERP: 405 code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n"; 406 code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" 407 << i << ", r_t_" << i << ");\n"; 408 break; 409 case CEED_EVAL_GRAD: 410 if (use_collograd_parallelization) { 411 code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*Q_1d];\n"; 412 code << " Interp" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_in_" << i << ",P_in_" << i << ",Q_1d>(data, r_u_" << i 413 << ", s_B_in_" << i << ", r_t_" << i << ");\n"; 414 } else { 415 CeedInt P_1d; 416 417 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 418 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 419 code << " CeedScalar r_t_" << i << "[num_comp_in_" << i << "*dim*Q_1d];\n"; 420 code << " Grad" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_in_" << i 421 << ",P_in_" << i << ",Q_1d>(data, r_u_" << i << ", s_B_in_" << i << ", s_G_in_" << i << ", r_t_" << i << ");\n"; 422 } 423 break; 424 case CEED_EVAL_WEIGHT: 425 code << " CeedScalar r_t_" << i << "[Q_1d];\n"; 426 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 427 CeedCallBackend(CeedBasisGetData(basis, &basis_data)); 428 data->W = basis_data->d_q_weight_1d; 429 code << " Weight" << (dim > 1 ? "Tensor" : "") << dim << "d<Q_1d>(data, W, r_t_" << i << ");\n"; 430 break; // No action 431 case CEED_EVAL_DIV: 432 break; // TODO: Not implemented 433 case CEED_EVAL_CURL: 434 break; // TODO: Not implemented 435 } 436 } 437 438 // TODO: put in a function + separate collograd logic 439 // Q function 440 code << "\n // -- Output field setup --\n"; 441 for (CeedInt i = 0; i < num_output_fields; i++) { 442 code << "\n // ---- Output field " << i << " ----\n"; 443 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 444 if (eval_mode == CEED_EVAL_GRAD) { 445 if (use_collograd_parallelization) { 446 // Accumulator for gradient slices 447 code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n"; 448 code << " for (CeedInt i = 0; i < num_comp_out_" << i << "; i++) {\n"; 449 code << " for (CeedInt j = 0; j < Q_1d; ++j) {\n"; 450 code << " r_tt_" << i << "[j + i*Q_1d] = 0.0;\n"; 451 code << " }\n"; 452 code << " }\n"; 453 } else { 454 code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*dim*Q_1d];\n"; 455 } 456 } 457 if (eval_mode == CEED_EVAL_NONE || eval_mode == CEED_EVAL_INTERP) { 458 code << " CeedScalar r_tt_" << i << "[num_comp_out_" << i << "*Q_1d];\n"; 459 } 460 } 461 // We treat quadrature points per slice in 3d to save registers 462 if (use_collograd_parallelization) { 463 code << "\n // Note: Using planes of 3D elements\n"; 464 code << "#pragma unroll\n"; 465 code << " for (CeedInt q = 0; q < Q_1d; q++) {\n"; 466 code << " // -- Input fields --\n"; 467 for (CeedInt i = 0; i < num_input_fields; i++) { 468 code << " // ---- Input field " << i << " ----\n"; 469 // Get elem_size, eval_mode, num_comp 470 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 471 // Basis action 472 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 473 switch (eval_mode) { 474 case CEED_EVAL_NONE: 475 bool is_strided; 476 477 code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n"; 478 479 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 480 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 481 if (!is_strided) { 482 CeedInt comp_stride; 483 484 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 485 code << " const CeedInt l_size_in_" << i << " = " << l_size << ";\n"; 486 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 487 code << " // CompStride: " << comp_stride << "\n"; 488 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 489 data->indices.inputs[i] = (CeedInt *)rstr_data->d_offsets; 490 code << " readSliceQuadsOffset" 491 << "3d<num_comp_in_" << i << ", " << comp_stride << ", Q_1d>(data, l_size_in_" << i << ", elem, q, indices.inputs[" << i << "], d_u_" 492 << i << ", r_q_" << i << ");\n"; 493 } else { 494 bool has_backend_strides; 495 CeedInt num_elem; 496 497 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 498 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 499 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 500 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 501 502 if (!has_backend_strides) { 503 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 504 } 505 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 506 code << " readSliceQuadsStrided" 507 << "3d<num_comp_in_" << i 508 << ",Q_1d" 509 "," 510 << strides[0] << "," << strides[1] << "," << strides[2] << ">(data, elem, q, d_u_" << i << ", r_q_" << i << ");\n"; 511 } 512 break; 513 case CEED_EVAL_INTERP: 514 code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "];\n"; 515 code << " for (CeedInt j = 0; j < num_comp_in_" << i << " ; ++j) {\n"; 516 code << " r_q_" << i << "[j] = r_t_" << i << "[q + j*Q_1d];\n"; 517 code << " }\n"; 518 break; 519 case CEED_EVAL_GRAD: 520 code << " CeedScalar r_q_" << i << "[num_comp_in_" << i << "*dim];\n"; 521 code << " gradCollo3d<num_comp_in_" << i << ",Q_1d>(data, q, r_t_" << i << ", s_G_in_" << i << ", r_q_" << i << ");\n"; 522 break; 523 case CEED_EVAL_WEIGHT: 524 code << " CeedScalar r_q_" << i << "[1];\n"; 525 code << " r_q_" << i << "[0] = r_t_" << i << "[q];\n"; 526 break; // No action 527 case CEED_EVAL_DIV: 528 break; // TODO: Not implemented 529 case CEED_EVAL_CURL: 530 break; // TODO: Not implemented 531 } 532 } 533 code << "\n // -- Output fields --\n"; 534 for (CeedInt i = 0; i < num_output_fields; i++) { 535 code << " // ---- Output field " << i << " ----\n"; 536 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 537 // Basis action 538 switch (eval_mode) { 539 case CEED_EVAL_NONE: 540 code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n"; 541 break; // No action 542 case CEED_EVAL_INTERP: 543 code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "];\n"; 544 break; 545 case CEED_EVAL_GRAD: 546 code << " CeedScalar r_qq_" << i << "[num_comp_out_" << i << "*dim];\n"; 547 break; 548 case CEED_EVAL_WEIGHT: 549 break; // Should not occur 550 case CEED_EVAL_DIV: 551 break; // TODO: Not implemented 552 case CEED_EVAL_CURL: 553 break; // TODO: Not implemented 554 } 555 } 556 } else { 557 code << "\n // Note: Using full elements\n"; 558 code << " // -- Input fields --\n"; 559 for (CeedInt i = 0; i < num_input_fields; i++) { 560 code << " // ---- Input field " << i << " ----\n"; 561 code << " CeedScalar* r_q_" << i << " = r_t_" << i << ";\n"; 562 } 563 code << " // -- Output fields --\n"; 564 for (CeedInt i = 0; i < num_output_fields; i++) { 565 code << " // ---- Output field " << i << " ----\n"; 566 code << " CeedScalar* r_qq_" << i << " = r_tt_" << i << ";\n"; 567 } 568 } 569 code << "\n // -- QFunction Inputs and outputs --\n"; 570 code << " CeedScalar* in[" << num_input_fields << "];\n"; 571 for (CeedInt i = 0; i < num_input_fields; i++) { 572 code << " // ---- Input field " << i << " ----\n"; 573 code << " in[" << i << "] = r_q_" << i << ";\n"; 574 } 575 code << " CeedScalar* out[" << num_output_fields << "];\n"; 576 for (CeedInt i = 0; i < num_output_fields; i++) { 577 code << " // ---- Output field " << i << " ----\n"; 578 code << " out[" << i << "] = r_qq_" << i << ";\n"; 579 } 580 code << "\n // -- Apply QFunction --\n"; 581 code << " " << qfunction_name << "(ctx, "; 582 if (dim != 3 || use_collograd_parallelization) { 583 code << "1"; 584 } else { 585 code << "Q_1d"; 586 } 587 code << ", in, out);\n"; 588 if (use_collograd_parallelization) { 589 code << " // -- Output fields --\n"; 590 for (CeedInt i = 0; i < num_output_fields; i++) { 591 code << " // ---- Output field " << i << " ----\n"; 592 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 593 // Basis action 594 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 595 switch (eval_mode) { 596 case CEED_EVAL_NONE: 597 code << " for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n"; 598 code << " r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n"; 599 code << " }\n"; 600 break; // No action 601 case CEED_EVAL_INTERP: 602 code << " for (CeedInt j = 0; j < num_comp_out_" << i << " ; ++j) {\n"; 603 code << " r_tt_" << i << "[q + j*Q_1d] = r_qq_" << i << "[j];\n"; 604 code << " }\n"; 605 break; 606 case CEED_EVAL_GRAD: 607 code << " gradColloTranspose3d<num_comp_out_" << i << ",Q_1d>(data, q, r_qq_" << i << ", s_G_out_" << i << ", r_tt_" << i << ");\n"; 608 break; 609 case CEED_EVAL_WEIGHT: 610 break; // Should not occur 611 case CEED_EVAL_DIV: 612 break; // TODO: Not implemented 613 case CEED_EVAL_CURL: 614 break; // TODO: Not implemented 615 } 616 } 617 code << " }\n"; 618 } 619 620 // Output basis apply if needed 621 // Generate the correct eval mode code for each output 622 code << "\n // -- Output field basis action and restrictions --\n"; 623 for (CeedInt i = 0; i < num_output_fields; i++) { 624 code << " // ---- Output field " << i << " ----\n"; 625 // Get elem_size, eval_mode, num_comp 626 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 627 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 628 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 629 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 630 // TODO put in a function 631 // Basis action 632 code << " // EvalMode: " << CeedEvalModes[eval_mode] << "\n"; 633 switch (eval_mode) { 634 case CEED_EVAL_NONE: 635 code << " CeedScalar* r_v_" << i << " = r_tt_" << i << ";\n"; 636 break; // No action 637 case CEED_EVAL_INTERP: 638 code << " CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n"; 639 code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i 640 << ", s_B_out_" << i << ", r_v_" << i << ");\n"; 641 break; 642 case CEED_EVAL_GRAD: 643 code << " CeedScalar r_v_" << i << "[num_comp_out_" << i << "*P_out_" << i << "];\n"; 644 if (use_collograd_parallelization) { 645 code << " InterpTranspose" << (dim > 1 ? "Tensor" : "") << dim << "d<num_comp_out_" << i << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i 646 << ", s_B_out_" << i << ", r_v_" << i << ");\n"; 647 } else { 648 CeedInt P_1d; 649 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 650 CeedCallBackend(CeedBasisGetNumNodes1D(basis, &P_1d)); 651 code << " GradTranspose" << (dim > 1 ? "Tensor" : "") << (dim == 3 && Q_1d >= P_1d ? "Collocated" : "") << dim << "d<num_comp_out_" << i 652 << ",P_out_" << i << ",Q_1d>(data, r_tt_" << i << ", s_B_out_" << i << ", s_G_out_" << i << ", r_v_" << i << ");\n"; 653 } 654 break; 655 // LCOV_EXCL_START 656 case CEED_EVAL_WEIGHT: { 657 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 658 break; // Should not occur 659 } 660 case CEED_EVAL_DIV: 661 case CEED_EVAL_CURL: { 662 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "%s not supported", CeedEvalModes[eval_mode]); 663 break; // Should not occur 664 } 665 // LCOV_EXCL_STOP 666 } 667 // TODO put in a function 668 // Restriction 669 bool is_strided; 670 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 671 if (!is_strided) { 672 CeedInt comp_stride; 673 674 CeedCallBackend(CeedElemRestrictionGetLVectorSize(elem_rstr, &l_size)); 675 code << " const CeedInt l_size_out_" << i << " = " << l_size << ";\n"; 676 CeedCallBackend(CeedElemRestrictionGetCompStride(elem_rstr, &comp_stride)); 677 code << " // CompStride: " << comp_stride << "\n"; 678 CeedCallBackend(CeedElemRestrictionGetData(elem_rstr, &rstr_data)); 679 data->indices.outputs[i] = (CeedInt *)rstr_data->d_offsets; 680 code << " writeDofsOffset" << dim << "d<num_comp_out_" << i << ", " << comp_stride << ", P_out_" << i << ">(data, l_size_out_" << i 681 << ", elem, indices.outputs[" << i << "], r_v_" << i << ", d_v_" << i << ");\n"; 682 } else { 683 bool has_backend_strides; 684 CeedInt num_elem; 685 686 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &has_backend_strides)); 687 CeedCallBackend(CeedElemRestrictionGetNumElements(elem_rstr, &num_elem)); 688 CeedInt strides[3] = {1, elem_size * num_elem, elem_size}; 689 690 if (!has_backend_strides) { 691 CeedCallBackend(CeedElemRestrictionGetStrides(elem_rstr, strides)); 692 } 693 code << " // Strides: {" << strides[0] << ", " << strides[1] << ", " << strides[2] << "}\n"; 694 code << " writeDofsStrided" << dim << "d<num_comp_out_" << i << ",P_out_" << i << "," << strides[0] << "," << strides[1] << "," << strides[2] 695 << ">(data, elem, r_v_" << i << ", d_v_" << i << ");\n"; 696 } 697 } 698 699 code << " }\n"; 700 code << "}\n"; 701 code << "// -----------------------------------------------------------------------------\n\n"; 702 703 // View kernel for debugging 704 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "Generated Operator Kernels:\n"); 705 CeedDebug(ceed, code.str().c_str()); 706 707 CeedCallBackend(CeedCompile_Cuda(ceed, code.str().c_str(), &data->module, 1, "T_1D", CeedIntMax(Q_1d, data->max_P_1d))); 708 CeedCallBackend(CeedGetKernel_Cuda(ceed, data->module, operator_name.c_str(), &data->op)); 709 710 CeedCallBackend(CeedOperatorSetSetupDone(op)); 711 return CEED_ERROR_SUCCESS; 712 } 713 714 //------------------------------------------------------------------------------ 715