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 #include <ceed.h> 9 #include <ceed/backend.h> 10 #include <ceed/jit-tools.h> 11 #include <assert.h> 12 #include <cuda.h> 13 #include <cuda_runtime.h> 14 #include <stdbool.h> 15 #include <string.h> 16 17 #include "../cuda/ceed-cuda-common.h" 18 #include "../cuda/ceed-cuda-compile.h" 19 #include "ceed-cuda-ref.h" 20 21 //------------------------------------------------------------------------------ 22 // Destroy operator 23 //------------------------------------------------------------------------------ 24 static int CeedOperatorDestroy_Cuda(CeedOperator op) { 25 CeedOperator_Cuda *impl; 26 27 CeedCallBackend(CeedOperatorGetData(op, &impl)); 28 29 // Apply data 30 CeedCallBackend(CeedFree(&impl->num_points)); 31 CeedCallBackend(CeedFree(&impl->skip_rstr_in)); 32 CeedCallBackend(CeedFree(&impl->skip_rstr_out)); 33 CeedCallBackend(CeedFree(&impl->apply_add_basis_out)); 34 for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) { 35 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i])); 36 } 37 CeedCallBackend(CeedFree(&impl->e_vecs)); 38 CeedCallBackend(CeedFree(&impl->input_states)); 39 40 for (CeedInt i = 0; i < impl->num_inputs; i++) { 41 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i])); 42 } 43 CeedCallBackend(CeedFree(&impl->q_vecs_in)); 44 45 for (CeedInt i = 0; i < impl->num_outputs; i++) { 46 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 47 } 48 CeedCallBackend(CeedFree(&impl->q_vecs_out)); 49 CeedCallBackend(CeedVectorDestroy(&impl->point_coords_elem)); 50 51 // QFunction assembly data 52 for (CeedInt i = 0; i < impl->num_active_in; i++) { 53 CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i])); 54 } 55 CeedCallBackend(CeedFree(&impl->qf_active_in)); 56 57 // Diag data 58 if (impl->diag) { 59 Ceed ceed; 60 61 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 62 if (impl->diag->module) { 63 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module)); 64 } 65 if (impl->diag->module_point_block) { 66 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block)); 67 } 68 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in)); 69 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out)); 70 CeedCallCuda(ceed, cudaFree(impl->diag->d_identity)); 71 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in)); 72 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out)); 73 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in)); 74 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out)); 75 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in)); 76 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out)); 77 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in)); 78 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out)); 79 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr)); 80 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr)); 81 CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag)); 82 CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag)); 83 } 84 CeedCallBackend(CeedFree(&impl->diag)); 85 86 if (impl->asmb) { 87 Ceed ceed; 88 89 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 90 CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module)); 91 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in)); 92 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out)); 93 } 94 CeedCallBackend(CeedFree(&impl->asmb)); 95 96 CeedCallBackend(CeedFree(&impl)); 97 return CEED_ERROR_SUCCESS; 98 } 99 100 //------------------------------------------------------------------------------ 101 // Setup infields or outfields 102 //------------------------------------------------------------------------------ 103 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, bool is_at_points, bool *skip_rstr, bool *apply_add_basis, 104 CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e, CeedInt num_fields, CeedInt Q, CeedInt num_elem) { 105 Ceed ceed; 106 CeedQFunctionField *qf_fields; 107 CeedOperatorField *op_fields; 108 109 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 110 if (is_input) { 111 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 112 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 113 } else { 114 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 115 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 116 } 117 118 // Loop over fields 119 for (CeedInt i = 0; i < num_fields; i++) { 120 bool is_strided = false, skip_restriction = false; 121 CeedSize q_size; 122 CeedInt size; 123 CeedEvalMode eval_mode; 124 CeedBasis basis; 125 126 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 127 if (eval_mode != CEED_EVAL_WEIGHT) { 128 CeedElemRestriction elem_rstr; 129 130 // Check whether this field can skip the element restriction: 131 // Must be passive input, with eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND. 132 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr)); 133 134 // First, check whether the field is input or output: 135 if (is_input) { 136 CeedVector vec; 137 138 // Check for passive input 139 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 140 if (vec != CEED_VECTOR_ACTIVE) { 141 // Check eval_mode 142 if (eval_mode == CEED_EVAL_NONE) { 143 // Check for strided restriction 144 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 145 if (is_strided) { 146 // Check if vector is already in preferred backend ordering 147 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_restriction)); 148 } 149 } 150 } 151 } 152 if (skip_restriction) { 153 // We do not need an E-Vector, but will use the input field vector's data directly in the operator application. 154 e_vecs[i + start_e] = NULL; 155 } else { 156 CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i + start_e])); 157 } 158 } 159 160 switch (eval_mode) { 161 case CEED_EVAL_NONE: 162 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 163 q_size = (CeedSize)num_elem * Q * size; 164 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 165 break; 166 case CEED_EVAL_INTERP: 167 case CEED_EVAL_GRAD: 168 case CEED_EVAL_DIV: 169 case CEED_EVAL_CURL: 170 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 171 q_size = (CeedSize)num_elem * Q * size; 172 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 173 break; 174 case CEED_EVAL_WEIGHT: // Only on input fields 175 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 176 q_size = (CeedSize)num_elem * Q; 177 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 178 if (is_at_points) { 179 CeedInt num_points[num_elem]; 180 181 for (CeedInt i = 0; i < num_elem; i++) num_points[i] = Q; 182 CeedCallBackend( 183 CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, CEED_VECTOR_NONE, q_vecs[i])); 184 } else { 185 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i])); 186 } 187 break; 188 } 189 } 190 // Drop duplicate restrictions 191 if (is_input) { 192 for (CeedInt i = 0; i < num_fields; i++) { 193 CeedVector vec_i; 194 CeedElemRestriction rstr_i; 195 196 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 197 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 198 for (CeedInt j = i + 1; j < num_fields; j++) { 199 CeedVector vec_j; 200 CeedElemRestriction rstr_j; 201 202 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 203 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 204 if (vec_i == vec_j && rstr_i == rstr_j) { 205 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i + start_e], &e_vecs[j + start_e])); 206 skip_rstr[j] = true; 207 } 208 } 209 } 210 } else { 211 for (CeedInt i = num_fields - 1; i >= 0; i--) { 212 CeedVector vec_i; 213 CeedElemRestriction rstr_i; 214 215 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 216 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 217 for (CeedInt j = i - 1; j >= 0; j--) { 218 CeedVector vec_j; 219 CeedElemRestriction rstr_j; 220 221 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 222 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 223 if (vec_i == vec_j && rstr_i == rstr_j) { 224 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i + start_e], &e_vecs[j + start_e])); 225 skip_rstr[j] = true; 226 apply_add_basis[i] = true; 227 } 228 } 229 } 230 } 231 return CEED_ERROR_SUCCESS; 232 } 233 234 //------------------------------------------------------------------------------ 235 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 236 //------------------------------------------------------------------------------ 237 static int CeedOperatorSetup_Cuda(CeedOperator op) { 238 Ceed ceed; 239 bool is_setup_done; 240 CeedInt Q, num_elem, num_input_fields, num_output_fields; 241 CeedQFunctionField *qf_input_fields, *qf_output_fields; 242 CeedQFunction qf; 243 CeedOperatorField *op_input_fields, *op_output_fields; 244 CeedOperator_Cuda *impl; 245 246 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 247 if (is_setup_done) return CEED_ERROR_SUCCESS; 248 249 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 250 CeedCallBackend(CeedOperatorGetData(op, &impl)); 251 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 252 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 253 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 254 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 255 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 256 257 // Allocate 258 CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs)); 259 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_in)); 260 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_out)); 261 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->apply_add_basis_out)); 262 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->input_states)); 263 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in)); 264 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out)); 265 impl->num_inputs = num_input_fields; 266 impl->num_outputs = num_output_fields; 267 268 // Set up infield and outfield e-vecs and q-vecs 269 // Infields 270 CeedCallBackend( 271 CeedOperatorSetupFields_Cuda(qf, op, true, false, impl->skip_rstr_in, NULL, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem)); 272 // Outfields 273 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, false, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs, impl->q_vecs_out, 274 num_input_fields, num_output_fields, Q, num_elem)); 275 276 // Reuse active e-vecs where able 277 { 278 CeedInt num_used = 0; 279 CeedElemRestriction *rstr_used = NULL; 280 281 for (CeedInt i = 0; i < num_input_fields; i++) { 282 bool is_used = false; 283 CeedVector vec_i; 284 CeedElemRestriction rstr_i; 285 286 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 287 if (vec_i != CEED_VECTOR_ACTIVE) continue; 288 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 289 for (CeedInt j = 0; j < num_used; j++) { 290 if (rstr_i == rstr_used[i]) is_used = true; 291 } 292 if (is_used) continue; 293 num_used++; 294 if (num_used == 1) CeedCallBackend(CeedCalloc(num_used, &rstr_used)); 295 else CeedCallBackend(CeedRealloc(num_used, &rstr_used)); 296 rstr_used[num_used - 1] = rstr_i; 297 for (CeedInt j = num_output_fields - 1; j >= 0; j--) { 298 CeedEvalMode eval_mode; 299 CeedVector vec_j; 300 CeedElemRestriction rstr_j; 301 302 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j)); 303 if (vec_j != CEED_VECTOR_ACTIVE) continue; 304 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode)); 305 if (eval_mode == CEED_EVAL_NONE) continue; 306 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j)); 307 if (rstr_i == rstr_j) { 308 CeedCallBackend(CeedVectorReferenceCopy(impl->e_vecs[i], &impl->e_vecs[j + impl->num_inputs])); 309 } 310 } 311 } 312 CeedCallBackend(CeedFree(&rstr_used)); 313 } 314 impl->has_shared_e_vecs = true; 315 CeedCallBackend(CeedOperatorSetSetupDone(op)); 316 return CEED_ERROR_SUCCESS; 317 } 318 319 //------------------------------------------------------------------------------ 320 // Restrict Operator Inputs 321 //------------------------------------------------------------------------------ 322 static inline int CeedOperatorInputRestrict_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 323 CeedVector in_vec, const bool skip_active, CeedScalar **e_data, CeedOperator_Cuda *impl, 324 CeedRequest *request) { 325 CeedEvalMode eval_mode; 326 CeedVector vec; 327 CeedElemRestriction elem_rstr; 328 329 // Get input vector 330 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 331 if (vec == CEED_VECTOR_ACTIVE) { 332 if (skip_active) return CEED_ERROR_SUCCESS; 333 else vec = in_vec; 334 } 335 336 // Restriction action 337 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 338 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 339 } else { 340 // Get input vector 341 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 342 // Get input element restriction 343 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_field, &elem_rstr)); 344 if (vec == CEED_VECTOR_ACTIVE) vec = in_vec; 345 // Restrict, if necessary 346 if (!impl->e_vecs[input_field]) { 347 // No restriction for this field; read data directly from vec. 348 CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)e_data)); 349 } else { 350 uint64_t state; 351 352 CeedCallBackend(CeedVectorGetState(vec, &state)); 353 if ((state != impl->input_states[input_field] || vec == in_vec) && !impl->skip_rstr_in[input_field]) { 354 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[input_field], request)); 355 } 356 impl->input_states[input_field] = state; 357 // Get evec 358 CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[input_field], CEED_MEM_DEVICE, (const CeedScalar **)e_data)); 359 } 360 } 361 return CEED_ERROR_SUCCESS; 362 } 363 364 //------------------------------------------------------------------------------ 365 // Input Basis Action 366 //------------------------------------------------------------------------------ 367 static inline int CeedOperatorInputBasis_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 368 CeedInt num_elem, const bool skip_active, CeedScalar *e_data, CeedOperator_Cuda *impl) { 369 CeedEvalMode eval_mode; 370 371 // Skip active input 372 if (skip_active) { 373 CeedVector vec; 374 375 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 376 if (vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 377 } 378 379 // Basis action 380 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 381 switch (eval_mode) { 382 case CEED_EVAL_NONE: 383 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[input_field], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data)); 384 break; 385 case CEED_EVAL_INTERP: 386 case CEED_EVAL_GRAD: 387 case CEED_EVAL_DIV: 388 case CEED_EVAL_CURL: { 389 CeedBasis basis; 390 391 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis)); 392 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs[input_field], impl->q_vecs_in[input_field])); 393 break; 394 } 395 case CEED_EVAL_WEIGHT: 396 break; // No action 397 } 398 return CEED_ERROR_SUCCESS; 399 } 400 401 //------------------------------------------------------------------------------ 402 // Restore Input Vectors 403 //------------------------------------------------------------------------------ 404 static inline int CeedOperatorInputRestore_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 405 const bool skip_active, CeedScalar **e_data, CeedOperator_Cuda *impl) { 406 CeedEvalMode eval_mode; 407 CeedVector vec; 408 409 // Skip active input 410 if (skip_active) { 411 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 412 if (vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 413 } 414 415 // Restore e-vec 416 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 417 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 418 } else { 419 if (!impl->e_vecs[input_field]) { // This was a skip_restriction case 420 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 421 CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)e_data)); 422 } else { 423 CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[input_field], (const CeedScalar **)e_data)); 424 } 425 } 426 return CEED_ERROR_SUCCESS; 427 } 428 429 //------------------------------------------------------------------------------ 430 // Apply and add to output 431 //------------------------------------------------------------------------------ 432 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 433 CeedInt Q, num_elem, elem_size, num_input_fields, num_output_fields, size; 434 CeedScalar *e_data[2 * CEED_FIELD_MAX] = {NULL}; 435 CeedQFunctionField *qf_input_fields, *qf_output_fields; 436 CeedQFunction qf; 437 CeedOperatorField *op_input_fields, *op_output_fields; 438 CeedOperator_Cuda *impl; 439 440 CeedCallBackend(CeedOperatorGetData(op, &impl)); 441 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 442 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 443 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 444 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 445 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 446 447 // Setup 448 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 449 450 // Process inputs 451 for (CeedInt i = 0; i < num_input_fields; i++) { 452 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, false, &e_data[i], impl, request)); 453 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, false, e_data[i], impl)); 454 } 455 456 // Output pointers, as necessary 457 for (CeedInt i = 0; i < num_output_fields; i++) { 458 CeedEvalMode eval_mode; 459 460 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 461 if (eval_mode == CEED_EVAL_NONE) { 462 // Set the output Q-Vector to use the E-Vector data directly. 463 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields])); 464 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields])); 465 } 466 } 467 468 // Q function 469 CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out)); 470 471 // Restore input arrays 472 for (CeedInt i = 0; i < num_input_fields; i++) { 473 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, false, &e_data[i], impl)); 474 } 475 476 // Output basis apply if needed 477 for (CeedInt i = 0; i < num_output_fields; i++) { 478 CeedEvalMode eval_mode; 479 CeedElemRestriction elem_rstr; 480 CeedBasis basis; 481 482 // Get elem_size, eval_mode, size 483 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 484 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 485 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 486 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 487 // Basis action 488 switch (eval_mode) { 489 case CEED_EVAL_NONE: 490 break; // No action 491 case CEED_EVAL_INTERP: 492 case CEED_EVAL_GRAD: 493 case CEED_EVAL_DIV: 494 case CEED_EVAL_CURL: 495 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 496 if (impl->apply_add_basis_out[i]) { 497 CeedCallBackend(CeedBasisApplyAdd(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs])); 498 } else { 499 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs])); 500 } 501 break; 502 // LCOV_EXCL_START 503 case CEED_EVAL_WEIGHT: { 504 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 505 // LCOV_EXCL_STOP 506 } 507 } 508 } 509 510 // Output restriction 511 for (CeedInt i = 0; i < num_output_fields; i++) { 512 CeedEvalMode eval_mode; 513 CeedVector vec; 514 CeedElemRestriction elem_rstr; 515 516 // Restore evec 517 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 518 if (eval_mode == CEED_EVAL_NONE) { 519 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields])); 520 } 521 if (impl->skip_rstr_out[i]) continue; 522 // Get output vector 523 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 524 // Restrict 525 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 526 // Active 527 if (vec == CEED_VECTOR_ACTIVE) vec = out_vec; 528 529 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request)); 530 } 531 return CEED_ERROR_SUCCESS; 532 } 533 534 //------------------------------------------------------------------------------ 535 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 536 //------------------------------------------------------------------------------ 537 static int CeedOperatorSetupAtPoints_Cuda(CeedOperator op) { 538 Ceed ceed; 539 bool is_setup_done; 540 CeedInt max_num_points = -1, num_elem, num_input_fields, num_output_fields; 541 CeedQFunctionField *qf_input_fields, *qf_output_fields; 542 CeedQFunction qf; 543 CeedOperatorField *op_input_fields, *op_output_fields; 544 CeedOperator_Cuda *impl; 545 546 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 547 if (is_setup_done) return CEED_ERROR_SUCCESS; 548 549 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 550 CeedCallBackend(CeedOperatorGetData(op, &impl)); 551 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 552 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 553 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 554 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 555 { 556 CeedElemRestriction rstr_points = NULL; 557 558 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, NULL)); 559 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(rstr_points, &max_num_points)); 560 CeedCallBackend(CeedCalloc(num_elem, &impl->num_points)); 561 for (CeedInt e = 0; e < num_elem; e++) { 562 CeedInt num_points_elem; 563 564 CeedCallBackend(CeedElemRestrictionGetNumPointsInElement(rstr_points, e, &num_points_elem)); 565 impl->num_points[e] = num_points_elem; 566 } 567 } 568 impl->max_num_points = max_num_points; 569 570 // Allocate 571 CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs)); 572 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_in)); 573 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_out)); 574 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->apply_add_basis_out)); 575 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->input_states)); 576 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in)); 577 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out)); 578 impl->num_inputs = num_input_fields; 579 impl->num_outputs = num_output_fields; 580 581 // Set up infield and outfield e-vecs and q-vecs 582 // Infields 583 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, true, impl->skip_rstr_in, NULL, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, 584 max_num_points, num_elem)); 585 // Outfields 586 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs, impl->q_vecs_out, 587 num_input_fields, num_output_fields, max_num_points, num_elem)); 588 589 // Reuse active e-vecs where able 590 { 591 CeedInt num_used = 0; 592 CeedElemRestriction *rstr_used = NULL; 593 594 for (CeedInt i = 0; i < num_input_fields; i++) { 595 bool is_used = false; 596 CeedVector vec_i; 597 CeedElemRestriction rstr_i; 598 599 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec_i)); 600 if (vec_i != CEED_VECTOR_ACTIVE) continue; 601 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &rstr_i)); 602 for (CeedInt j = 0; j < num_used; j++) { 603 if (rstr_i == rstr_used[i]) is_used = true; 604 } 605 if (is_used) continue; 606 num_used++; 607 if (num_used == 1) CeedCallBackend(CeedCalloc(num_used, &rstr_used)); 608 else CeedCallBackend(CeedRealloc(num_used, &rstr_used)); 609 rstr_used[num_used - 1] = rstr_i; 610 for (CeedInt j = num_output_fields - 1; j >= 0; j--) { 611 CeedEvalMode eval_mode; 612 CeedVector vec_j; 613 CeedElemRestriction rstr_j; 614 615 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec_j)); 616 if (vec_j != CEED_VECTOR_ACTIVE) continue; 617 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode)); 618 if (eval_mode == CEED_EVAL_NONE) continue; 619 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &rstr_j)); 620 if (rstr_i == rstr_j) { 621 CeedCallBackend(CeedVectorReferenceCopy(impl->e_vecs[i], &impl->e_vecs[j + impl->num_inputs])); 622 } 623 } 624 } 625 CeedCallBackend(CeedFree(&rstr_used)); 626 } 627 impl->has_shared_e_vecs = true; 628 CeedCallBackend(CeedOperatorSetSetupDone(op)); 629 return CEED_ERROR_SUCCESS; 630 } 631 632 //------------------------------------------------------------------------------ 633 // Input Basis Action AtPoints 634 //------------------------------------------------------------------------------ 635 static inline int CeedOperatorInputBasisAtPoints_Cuda(CeedOperatorField op_input_field, CeedQFunctionField qf_input_field, CeedInt input_field, 636 CeedInt num_elem, const CeedInt *num_points, const bool skip_active, CeedScalar *e_data, 637 CeedOperator_Cuda *impl) { 638 CeedEvalMode eval_mode; 639 640 // Skip active input 641 if (skip_active) { 642 CeedVector vec; 643 644 CeedCallBackend(CeedOperatorFieldGetVector(op_input_field, &vec)); 645 if (vec == CEED_VECTOR_ACTIVE) return CEED_ERROR_SUCCESS; 646 } 647 648 // Basis action 649 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_field, &eval_mode)); 650 switch (eval_mode) { 651 case CEED_EVAL_NONE: 652 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[input_field], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data)); 653 break; 654 case CEED_EVAL_INTERP: 655 case CEED_EVAL_GRAD: 656 case CEED_EVAL_DIV: 657 case CEED_EVAL_CURL: { 658 CeedBasis basis; 659 660 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_field, &basis)); 661 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, 662 impl->e_vecs[input_field], impl->q_vecs_in[input_field])); 663 break; 664 } 665 case CEED_EVAL_WEIGHT: 666 break; // No action 667 } 668 return CEED_ERROR_SUCCESS; 669 } 670 671 //------------------------------------------------------------------------------ 672 // Apply and add to output AtPoints 673 //------------------------------------------------------------------------------ 674 static int CeedOperatorApplyAddAtPoints_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 675 CeedInt max_num_points, *num_points, num_elem, elem_size, num_input_fields, num_output_fields, size; 676 CeedScalar *e_data[2 * CEED_FIELD_MAX] = {NULL}; 677 CeedQFunctionField *qf_input_fields, *qf_output_fields; 678 CeedQFunction qf; 679 CeedOperatorField *op_input_fields, *op_output_fields; 680 CeedOperator_Cuda *impl; 681 682 CeedCallBackend(CeedOperatorGetData(op, &impl)); 683 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 684 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 685 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 686 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 687 688 // Setup 689 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op)); 690 num_points = impl->num_points; 691 max_num_points = impl->max_num_points; 692 693 // Get point coordinates 694 if (!impl->point_coords_elem) { 695 CeedVector point_coords = NULL; 696 CeedElemRestriction rstr_points = NULL; 697 698 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords)); 699 CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem)); 700 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request)); 701 } 702 703 // Process inputs 704 for (CeedInt i = 0; i < num_input_fields; i++) { 705 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, in_vec, false, &e_data[i], impl, request)); 706 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, num_points, false, e_data[i], impl)); 707 } 708 709 // Output pointers, as necessary 710 for (CeedInt i = 0; i < num_output_fields; i++) { 711 CeedEvalMode eval_mode; 712 713 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 714 if (eval_mode == CEED_EVAL_NONE) { 715 // Set the output Q-Vector to use the E-Vector data directly. 716 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields])); 717 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields])); 718 } 719 } 720 721 // Q function 722 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out)); 723 724 // Restore input arrays 725 for (CeedInt i = 0; i < num_input_fields; i++) { 726 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, false, &e_data[i], impl)); 727 } 728 729 // Output basis apply if needed 730 for (CeedInt i = 0; i < num_output_fields; i++) { 731 CeedEvalMode eval_mode; 732 CeedElemRestriction elem_rstr; 733 CeedBasis basis; 734 735 // Get elem_size, eval_mode, size 736 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 737 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 738 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 739 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 740 // Basis action 741 switch (eval_mode) { 742 case CEED_EVAL_NONE: 743 break; // No action 744 case CEED_EVAL_INTERP: 745 case CEED_EVAL_GRAD: 746 case CEED_EVAL_DIV: 747 case CEED_EVAL_CURL: 748 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 749 if (impl->apply_add_basis_out[i]) { 750 CeedCallBackend(CeedBasisApplyAddAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, 751 impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs])); 752 } else { 753 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, impl->q_vecs_out[i], 754 impl->e_vecs[i + impl->num_inputs])); 755 } 756 break; 757 // LCOV_EXCL_START 758 case CEED_EVAL_WEIGHT: { 759 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 760 // LCOV_EXCL_STOP 761 } 762 } 763 } 764 765 // Output restriction 766 for (CeedInt i = 0; i < num_output_fields; i++) { 767 CeedEvalMode eval_mode; 768 CeedVector vec; 769 CeedElemRestriction elem_rstr; 770 771 // Restore evec 772 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 773 if (eval_mode == CEED_EVAL_NONE) { 774 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields])); 775 } 776 if (impl->skip_rstr_out[i]) continue; 777 // Get output vector 778 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 779 // Restrict 780 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 781 // Active 782 if (vec == CEED_VECTOR_ACTIVE) vec = out_vec; 783 784 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request)); 785 } 786 return CEED_ERROR_SUCCESS; 787 } 788 789 //------------------------------------------------------------------------------ 790 // Linear QFunction Assembly Core 791 //------------------------------------------------------------------------------ 792 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 793 CeedRequest *request) { 794 Ceed ceed, ceed_parent; 795 CeedInt num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size; 796 CeedScalar *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL}; 797 CeedVector *active_inputs; 798 CeedQFunctionField *qf_input_fields, *qf_output_fields; 799 CeedQFunction qf; 800 CeedOperatorField *op_input_fields, *op_output_fields; 801 CeedOperator_Cuda *impl; 802 803 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 804 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent)); 805 CeedCallBackend(CeedOperatorGetData(op, &impl)); 806 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 807 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 808 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 809 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 810 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 811 active_inputs = impl->qf_active_in; 812 num_active_in = impl->num_active_in, num_active_out = impl->num_active_out; 813 814 // Setup 815 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 816 817 // Input Evecs and Restriction 818 for (CeedInt i = 0; i < num_input_fields; i++) { 819 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, true, &e_data[i], impl, request)); 820 } 821 822 // Count number of active input fields 823 if (!num_active_in) { 824 for (CeedInt i = 0; i < num_input_fields; i++) { 825 CeedScalar *q_vec_array; 826 CeedVector vec; 827 828 // Get input vector 829 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 830 // Check if active input 831 if (vec == CEED_VECTOR_ACTIVE) { 832 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 833 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 834 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array)); 835 CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs)); 836 for (CeedInt field = 0; field < size; field++) { 837 CeedSize q_size = (CeedSize)Q * num_elem; 838 839 CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field])); 840 CeedCallBackend( 841 CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem])); 842 } 843 num_active_in += size; 844 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array)); 845 } 846 } 847 impl->num_active_in = num_active_in; 848 impl->qf_active_in = active_inputs; 849 } 850 851 // Count number of active output fields 852 if (!num_active_out) { 853 for (CeedInt i = 0; i < num_output_fields; i++) { 854 CeedVector vec; 855 856 // Get output vector 857 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 858 // Check if active output 859 if (vec == CEED_VECTOR_ACTIVE) { 860 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 861 num_active_out += size; 862 } 863 } 864 impl->num_active_out = num_active_out; 865 } 866 867 // Check sizes 868 CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 869 870 // Build objects if needed 871 if (build_objects) { 872 CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out; 873 CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */ 874 875 // Create output restriction 876 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out, 877 (CeedSize)num_active_in * (CeedSize)num_active_out * (CeedSize)num_elem * (CeedSize)Q, strides, 878 rstr)); 879 // Create assembled vector 880 CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled)); 881 } 882 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0)); 883 CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array)); 884 885 // Input basis apply 886 for (CeedInt i = 0; i < num_input_fields; i++) { 887 CeedCallBackend(CeedOperatorInputBasis_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, true, e_data[i], impl)); 888 } 889 890 // Assemble QFunction 891 for (CeedInt in = 0; in < num_active_in; in++) { 892 // Set Inputs 893 CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0)); 894 if (num_active_in > 1) { 895 CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0)); 896 } 897 // Set Outputs 898 for (CeedInt out = 0; out < num_output_fields; out++) { 899 CeedVector vec; 900 901 // Get output vector 902 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 903 // Check if active output 904 if (vec == CEED_VECTOR_ACTIVE) { 905 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array)); 906 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size)); 907 assembled_array += size * Q * num_elem; // Advance the pointer by the size of the output 908 } 909 } 910 // Apply QFunction 911 CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out)); 912 } 913 914 // Un-set output q-vecs to prevent accidental overwrite of Assembled 915 for (CeedInt out = 0; out < num_output_fields; out++) { 916 CeedVector vec; 917 918 // Get output vector 919 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 920 // Check if active output 921 if (vec == CEED_VECTOR_ACTIVE) { 922 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL)); 923 } 924 } 925 926 // Restore input arrays 927 for (CeedInt i = 0; i < num_input_fields; i++) { 928 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, true, &e_data[i], impl)); 929 } 930 931 // Restore output 932 CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array)); 933 return CEED_ERROR_SUCCESS; 934 } 935 936 //------------------------------------------------------------------------------ 937 // Assemble Linear QFunction 938 //------------------------------------------------------------------------------ 939 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 940 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request); 941 } 942 943 //------------------------------------------------------------------------------ 944 // Update Assembled Linear QFunction 945 //------------------------------------------------------------------------------ 946 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 947 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request); 948 } 949 950 //------------------------------------------------------------------------------ 951 // Assemble Diagonal Setup 952 //------------------------------------------------------------------------------ 953 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) { 954 Ceed ceed; 955 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 956 CeedInt q_comp, num_nodes, num_qpts; 957 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 958 CeedBasis basis_in = NULL, basis_out = NULL; 959 CeedQFunctionField *qf_fields; 960 CeedQFunction qf; 961 CeedOperatorField *op_fields; 962 CeedOperator_Cuda *impl; 963 964 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 965 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 966 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 967 968 // Determine active input basis 969 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 970 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 971 for (CeedInt i = 0; i < num_input_fields; i++) { 972 CeedVector vec; 973 974 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 975 if (vec == CEED_VECTOR_ACTIVE) { 976 CeedBasis basis; 977 CeedEvalMode eval_mode; 978 979 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 980 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, 981 "Backend does not implement operator diagonal assembly with multiple active bases"); 982 basis_in = basis; 983 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 984 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 985 if (eval_mode != CEED_EVAL_WEIGHT) { 986 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 987 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 988 for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode; 989 num_eval_modes_in += q_comp; 990 } 991 } 992 } 993 994 // Determine active output basis 995 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 996 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 997 for (CeedInt i = 0; i < num_output_fields; i++) { 998 CeedVector vec; 999 1000 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1001 if (vec == CEED_VECTOR_ACTIVE) { 1002 CeedBasis basis; 1003 CeedEvalMode eval_mode; 1004 1005 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 1006 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 1007 "Backend does not implement operator diagonal assembly with multiple active bases"); 1008 basis_out = basis; 1009 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1010 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1011 if (eval_mode != CEED_EVAL_WEIGHT) { 1012 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 1013 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 1014 for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode; 1015 num_eval_modes_out += q_comp; 1016 } 1017 } 1018 } 1019 1020 // Operator data struct 1021 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1022 CeedCallBackend(CeedCalloc(1, &impl->diag)); 1023 CeedOperatorDiag_Cuda *diag = impl->diag; 1024 1025 // Basis matrices 1026 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 1027 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 1028 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 1029 const CeedInt interp_bytes = num_nodes * num_qpts * sizeof(CeedScalar); 1030 const CeedInt eval_modes_bytes = sizeof(CeedEvalMode); 1031 bool has_eval_none = false; 1032 1033 // CEED_EVAL_NONE 1034 for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 1035 for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 1036 if (has_eval_none) { 1037 CeedScalar *identity = NULL; 1038 1039 CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity)); 1040 for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0; 1041 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes)); 1042 CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice)); 1043 CeedCallBackend(CeedFree(&identity)); 1044 } 1045 1046 // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL 1047 for (CeedInt in = 0; in < 2; in++) { 1048 CeedFESpace fespace; 1049 CeedBasis basis = in ? basis_in : basis_out; 1050 1051 CeedCallBackend(CeedBasisGetFESpace(basis, &fespace)); 1052 switch (fespace) { 1053 case CEED_FE_SPACE_H1: { 1054 CeedInt q_comp_interp, q_comp_grad; 1055 const CeedScalar *interp, *grad; 1056 CeedScalar *d_interp, *d_grad; 1057 1058 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1059 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 1060 1061 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1062 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1063 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1064 CeedCallBackend(CeedBasisGetGrad(basis, &grad)); 1065 CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad)); 1066 CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice)); 1067 if (in) { 1068 diag->d_interp_in = d_interp; 1069 diag->d_grad_in = d_grad; 1070 } else { 1071 diag->d_interp_out = d_interp; 1072 diag->d_grad_out = d_grad; 1073 } 1074 } break; 1075 case CEED_FE_SPACE_HDIV: { 1076 CeedInt q_comp_interp, q_comp_div; 1077 const CeedScalar *interp, *div; 1078 CeedScalar *d_interp, *d_div; 1079 1080 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1081 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 1082 1083 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1084 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1085 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1086 CeedCallBackend(CeedBasisGetDiv(basis, &div)); 1087 CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div)); 1088 CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice)); 1089 if (in) { 1090 diag->d_interp_in = d_interp; 1091 diag->d_div_in = d_div; 1092 } else { 1093 diag->d_interp_out = d_interp; 1094 diag->d_div_out = d_div; 1095 } 1096 } break; 1097 case CEED_FE_SPACE_HCURL: { 1098 CeedInt q_comp_interp, q_comp_curl; 1099 const CeedScalar *interp, *curl; 1100 CeedScalar *d_interp, *d_curl; 1101 1102 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 1103 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 1104 1105 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 1106 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 1107 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 1108 CeedCallBackend(CeedBasisGetCurl(basis, &curl)); 1109 CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl)); 1110 CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice)); 1111 if (in) { 1112 diag->d_interp_in = d_interp; 1113 diag->d_curl_in = d_curl; 1114 } else { 1115 diag->d_interp_out = d_interp; 1116 diag->d_curl_out = d_curl; 1117 } 1118 } break; 1119 } 1120 } 1121 1122 // Arrays of eval_modes 1123 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes)); 1124 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice)); 1125 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes)); 1126 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice)); 1127 CeedCallBackend(CeedFree(&eval_modes_in)); 1128 CeedCallBackend(CeedFree(&eval_modes_out)); 1129 return CEED_ERROR_SUCCESS; 1130 } 1131 1132 //------------------------------------------------------------------------------ 1133 // Assemble Diagonal Setup (Compilation) 1134 //------------------------------------------------------------------------------ 1135 static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) { 1136 Ceed ceed; 1137 char *diagonal_kernel_source; 1138 const char *diagonal_kernel_path; 1139 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 1140 CeedInt num_comp, q_comp, num_nodes, num_qpts; 1141 CeedBasis basis_in = NULL, basis_out = NULL; 1142 CeedQFunctionField *qf_fields; 1143 CeedQFunction qf; 1144 CeedOperatorField *op_fields; 1145 CeedOperator_Cuda *impl; 1146 1147 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1148 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1149 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 1150 1151 // Determine active input basis 1152 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 1153 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 1154 for (CeedInt i = 0; i < num_input_fields; i++) { 1155 CeedVector vec; 1156 1157 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1158 if (vec == CEED_VECTOR_ACTIVE) { 1159 CeedEvalMode eval_mode; 1160 1161 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in)); 1162 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1163 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 1164 if (eval_mode != CEED_EVAL_WEIGHT) { 1165 num_eval_modes_in += q_comp; 1166 } 1167 } 1168 } 1169 1170 // Determine active output basis 1171 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 1172 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1173 for (CeedInt i = 0; i < num_output_fields; i++) { 1174 CeedVector vec; 1175 1176 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 1177 if (vec == CEED_VECTOR_ACTIVE) { 1178 CeedEvalMode eval_mode; 1179 1180 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out)); 1181 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1182 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1183 if (eval_mode != CEED_EVAL_WEIGHT) { 1184 num_eval_modes_out += q_comp; 1185 } 1186 } 1187 } 1188 1189 // Operator data struct 1190 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1191 CeedOperatorDiag_Cuda *diag = impl->diag; 1192 1193 // Assemble kernel 1194 CUmodule *module = is_point_block ? &diag->module_point_block : &diag->module; 1195 CeedInt elems_per_block = 1; 1196 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 1197 CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp)); 1198 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 1199 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 1200 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path)); 1201 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n"); 1202 CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source)); 1203 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n"); 1204 CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 1205 num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE", 1206 use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block)); 1207 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal)); 1208 CeedCallBackend(CeedFree(&diagonal_kernel_path)); 1209 CeedCallBackend(CeedFree(&diagonal_kernel_source)); 1210 return CEED_ERROR_SUCCESS; 1211 } 1212 1213 //------------------------------------------------------------------------------ 1214 // Assemble Diagonal Core 1215 //------------------------------------------------------------------------------ 1216 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) { 1217 Ceed ceed; 1218 CeedInt num_elem, num_nodes; 1219 CeedScalar *elem_diag_array; 1220 const CeedScalar *assembled_qf_array; 1221 CeedVector assembled_qf = NULL, elem_diag; 1222 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr; 1223 CeedOperator_Cuda *impl; 1224 1225 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1226 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1227 1228 // Assemble QFunction 1229 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request)); 1230 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 1231 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 1232 1233 // Setup 1234 if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op)); 1235 CeedOperatorDiag_Cuda *diag = impl->diag; 1236 1237 assert(diag != NULL); 1238 1239 // Assemble kernel if needed 1240 if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) { 1241 CeedSize assembled_length, assembled_qf_length; 1242 CeedInt use_ceedsize_idx = 0; 1243 CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length)); 1244 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1245 if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1246 1247 CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block)); 1248 } 1249 1250 // Restriction and diagonal vector 1251 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 1252 CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND, 1253 "Cannot assemble operator diagonal with different input and output active element restrictions"); 1254 if (!is_point_block && !diag->diag_rstr) { 1255 CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr)); 1256 CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag)); 1257 } else if (is_point_block && !diag->point_block_diag_rstr) { 1258 CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr)); 1259 CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag)); 1260 } 1261 diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr; 1262 elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag; 1263 CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0)); 1264 1265 // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers 1266 CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes)); 1267 if (num_nodes > 0) { 1268 // Assemble element operator diagonals 1269 CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem)); 1270 CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array)); 1271 1272 // Compute the diagonal of B^T D B 1273 CeedInt elems_per_block = 1; 1274 CeedInt grid = CeedDivUpInt(num_elem, elems_per_block); 1275 void *args[] = {(void *)&num_elem, &diag->d_identity, &diag->d_interp_in, &diag->d_grad_in, &diag->d_div_in, 1276 &diag->d_curl_in, &diag->d_interp_out, &diag->d_grad_out, &diag->d_div_out, &diag->d_curl_out, 1277 &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array}; 1278 1279 if (is_point_block) { 1280 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args)); 1281 } else { 1282 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args)); 1283 } 1284 1285 // Restore arrays 1286 CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array)); 1287 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 1288 } 1289 1290 // Assemble local operator diagonal 1291 CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request)); 1292 1293 // Cleanup 1294 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1295 return CEED_ERROR_SUCCESS; 1296 } 1297 1298 //------------------------------------------------------------------------------ 1299 // Assemble Linear Diagonal 1300 //------------------------------------------------------------------------------ 1301 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1302 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false)); 1303 return CEED_ERROR_SUCCESS; 1304 } 1305 1306 //------------------------------------------------------------------------------ 1307 // Assemble Linear Point Block Diagonal 1308 //------------------------------------------------------------------------------ 1309 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1310 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true)); 1311 return CEED_ERROR_SUCCESS; 1312 } 1313 1314 //------------------------------------------------------------------------------ 1315 // Single Operator Assembly Setup 1316 //------------------------------------------------------------------------------ 1317 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) { 1318 Ceed ceed; 1319 Ceed_Cuda *cuda_data; 1320 char *assembly_kernel_source; 1321 const char *assembly_kernel_path; 1322 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 1323 CeedInt elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp; 1324 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 1325 CeedElemRestriction rstr_in = NULL, rstr_out = NULL; 1326 CeedBasis basis_in = NULL, basis_out = NULL; 1327 CeedQFunctionField *qf_fields; 1328 CeedQFunction qf; 1329 CeedOperatorField *input_fields, *output_fields; 1330 CeedOperator_Cuda *impl; 1331 1332 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1333 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1334 1335 // Get intput and output fields 1336 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 1337 1338 // Determine active input basis eval mode 1339 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1340 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 1341 for (CeedInt i = 0; i < num_input_fields; i++) { 1342 CeedVector vec; 1343 1344 CeedCallBackend(CeedOperatorFieldGetVector(input_fields[i], &vec)); 1345 if (vec == CEED_VECTOR_ACTIVE) { 1346 CeedBasis basis; 1347 CeedEvalMode eval_mode; 1348 1349 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis)); 1350 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases"); 1351 basis_in = basis; 1352 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in)); 1353 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 1354 if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in; 1355 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in)); 1356 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1357 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 1358 if (eval_mode != CEED_EVAL_WEIGHT) { 1359 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 1360 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 1361 for (CeedInt d = 0; d < q_comp; d++) { 1362 eval_modes_in[num_eval_modes_in + d] = eval_mode; 1363 } 1364 num_eval_modes_in += q_comp; 1365 } 1366 } 1367 } 1368 1369 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 1370 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1371 for (CeedInt i = 0; i < num_output_fields; i++) { 1372 CeedVector vec; 1373 1374 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec)); 1375 if (vec == CEED_VECTOR_ACTIVE) { 1376 CeedBasis basis; 1377 CeedEvalMode eval_mode; 1378 1379 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis)); 1380 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 1381 "Backend does not implement operator assembly with multiple active bases"); 1382 basis_out = basis; 1383 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out)); 1384 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1385 if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out; 1386 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out)); 1387 CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED, 1388 "Active input and output bases must have the same number of quadrature points"); 1389 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1390 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1391 if (eval_mode != CEED_EVAL_WEIGHT) { 1392 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 1393 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 1394 for (CeedInt d = 0; d < q_comp; d++) { 1395 eval_modes_out[num_eval_modes_out + d] = eval_mode; 1396 } 1397 num_eval_modes_out += q_comp; 1398 } 1399 } 1400 } 1401 CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs"); 1402 1403 CeedCallBackend(CeedCalloc(1, &impl->asmb)); 1404 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1405 asmb->elems_per_block = 1; 1406 asmb->block_size_x = elem_size_in; 1407 asmb->block_size_y = elem_size_out; 1408 1409 CeedCallBackend(CeedGetData(ceed, &cuda_data)); 1410 bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock; 1411 1412 if (fallback) { 1413 // Use fallback kernel with 1D threadblock 1414 asmb->block_size_y = 1; 1415 } 1416 1417 // Compile kernels 1418 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in)); 1419 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out)); 1420 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path)); 1421 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n"); 1422 CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source)); 1423 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n"); 1424 CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 1425 num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in, 1426 "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE", 1427 asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y, 1428 "USE_CEEDSIZE", use_ceedsize_idx)); 1429 CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble)); 1430 CeedCallBackend(CeedFree(&assembly_kernel_path)); 1431 CeedCallBackend(CeedFree(&assembly_kernel_source)); 1432 1433 // Load into B_in, in order that they will be used in eval_modes_in 1434 { 1435 const CeedInt in_bytes = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar); 1436 CeedInt d_in = 0; 1437 CeedEvalMode eval_modes_in_prev = CEED_EVAL_NONE; 1438 bool has_eval_none = false; 1439 CeedScalar *identity = NULL; 1440 1441 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1442 has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 1443 } 1444 if (has_eval_none) { 1445 CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity)); 1446 for (CeedInt i = 0; i < (elem_size_in < num_qpts_in ? elem_size_in : num_qpts_in); i++) identity[i * elem_size_in + i] = 1.0; 1447 } 1448 1449 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes)); 1450 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1451 const CeedScalar *h_B_in; 1452 1453 CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in)); 1454 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp)); 1455 if (q_comp > 1) { 1456 if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0; 1457 else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in]; 1458 } 1459 eval_modes_in_prev = eval_modes_in[i]; 1460 1461 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_in[i * elem_size_in * num_qpts_in], h_B_in, elem_size_in * num_qpts_in * sizeof(CeedScalar), 1462 cudaMemcpyHostToDevice)); 1463 } 1464 1465 if (identity) { 1466 CeedCallBackend(CeedFree(&identity)); 1467 } 1468 } 1469 1470 // Load into B_out, in order that they will be used in eval_modes_out 1471 { 1472 const CeedInt out_bytes = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar); 1473 CeedInt d_out = 0; 1474 CeedEvalMode eval_modes_out_prev = CEED_EVAL_NONE; 1475 bool has_eval_none = false; 1476 CeedScalar *identity = NULL; 1477 1478 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1479 has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 1480 } 1481 if (has_eval_none) { 1482 CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity)); 1483 for (CeedInt i = 0; i < (elem_size_out < num_qpts_out ? elem_size_out : num_qpts_out); i++) identity[i * elem_size_out + i] = 1.0; 1484 } 1485 1486 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes)); 1487 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1488 const CeedScalar *h_B_out; 1489 1490 CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out)); 1491 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp)); 1492 if (q_comp > 1) { 1493 if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0; 1494 else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out]; 1495 } 1496 eval_modes_out_prev = eval_modes_out[i]; 1497 1498 CeedCallCuda(ceed, cudaMemcpy(&asmb->d_B_out[i * elem_size_out * num_qpts_out], h_B_out, elem_size_out * num_qpts_out * sizeof(CeedScalar), 1499 cudaMemcpyHostToDevice)); 1500 } 1501 1502 if (identity) { 1503 CeedCallBackend(CeedFree(&identity)); 1504 } 1505 } 1506 return CEED_ERROR_SUCCESS; 1507 } 1508 1509 //------------------------------------------------------------------------------ 1510 // Assemble matrix data for COO matrix of assembled operator. 1511 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1512 // 1513 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval 1514 // modes). 1515 // TODO: allow multiple active input restrictions/basis objects 1516 //------------------------------------------------------------------------------ 1517 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) { 1518 Ceed ceed; 1519 CeedSize values_length = 0, assembled_qf_length = 0; 1520 CeedInt use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out; 1521 CeedScalar *values_array; 1522 const CeedScalar *assembled_qf_array; 1523 CeedVector assembled_qf = NULL; 1524 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out; 1525 CeedRestrictionType rstr_type_in, rstr_type_out; 1526 const bool *orients_in = NULL, *orients_out = NULL; 1527 const CeedInt8 *curl_orients_in = NULL, *curl_orients_out = NULL; 1528 CeedOperator_Cuda *impl; 1529 1530 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1531 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1532 1533 // Assemble QFunction 1534 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE)); 1535 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 1536 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 1537 1538 CeedCallBackend(CeedVectorGetLength(values, &values_length)); 1539 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1540 if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1541 1542 // Setup 1543 if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx)); 1544 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1545 1546 assert(asmb != NULL); 1547 1548 // Assemble element operator 1549 CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array)); 1550 values_array += offset; 1551 1552 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 1553 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in)); 1554 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 1555 1556 CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in)); 1557 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1558 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in)); 1559 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1560 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in)); 1561 } 1562 1563 if (rstr_in != rstr_out) { 1564 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out)); 1565 CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED, 1566 "Active input and output operator restrictions must have the same number of elements"); 1567 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1568 1569 CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out)); 1570 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1571 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out)); 1572 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1573 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out)); 1574 } 1575 } else { 1576 elem_size_out = elem_size_in; 1577 orients_out = orients_in; 1578 curl_orients_out = curl_orients_in; 1579 } 1580 1581 // Compute B^T D B 1582 CeedInt shared_mem = 1583 ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) * 1584 sizeof(CeedScalar); 1585 CeedInt grid = CeedDivUpInt(num_elem_in, asmb->elems_per_block); 1586 void *args[] = {(void *)&num_elem_in, &asmb->d_B_in, &asmb->d_B_out, &orients_in, &curl_orients_in, 1587 &orients_out, &curl_orients_out, &assembled_qf_array, &values_array}; 1588 1589 CeedCallBackend( 1590 CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args)); 1591 1592 // Restore arrays 1593 CeedCallBackend(CeedVectorRestoreArray(values, &values_array)); 1594 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 1595 1596 // Cleanup 1597 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1598 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1599 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in)); 1600 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1601 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in)); 1602 } 1603 if (rstr_in != rstr_out) { 1604 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1605 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out)); 1606 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1607 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out)); 1608 } 1609 } 1610 return CEED_ERROR_SUCCESS; 1611 } 1612 1613 //------------------------------------------------------------------------------ 1614 // Assemble Linear QFunction AtPoints 1615 //------------------------------------------------------------------------------ 1616 static int CeedOperatorLinearAssembleQFunctionAtPoints_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 1617 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "Backend does not implement CeedOperatorLinearAssembleQFunction"); 1618 } 1619 1620 //------------------------------------------------------------------------------ 1621 // Assemble Linear Diagonal AtPoints 1622 //------------------------------------------------------------------------------ 1623 static int CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 1624 CeedInt max_num_points, *num_points, num_elem, num_input_fields, num_output_fields; 1625 CeedScalar *e_data[2 * CEED_FIELD_MAX] = {NULL}; 1626 CeedQFunctionField *qf_input_fields, *qf_output_fields; 1627 CeedQFunction qf; 1628 CeedOperatorField *op_input_fields, *op_output_fields; 1629 CeedOperator_Cuda *impl; 1630 1631 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1632 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 1633 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 1634 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 1635 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 1636 1637 // Setup 1638 CeedCallBackend(CeedOperatorSetupAtPoints_Cuda(op)); 1639 num_points = impl->num_points; 1640 max_num_points = impl->max_num_points; 1641 1642 // Create separate output e-vecs 1643 if (impl->has_shared_e_vecs) { 1644 for (CeedInt i = 0; i < impl->num_outputs; i++) { 1645 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 1646 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[impl->num_inputs + i])); 1647 } 1648 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, true, impl->skip_rstr_out, impl->apply_add_basis_out, impl->e_vecs, impl->q_vecs_out, 1649 num_input_fields, num_output_fields, max_num_points, num_elem)); 1650 } 1651 impl->has_shared_e_vecs = false; 1652 1653 // Get point coordinates 1654 if (!impl->point_coords_elem) { 1655 CeedVector point_coords = NULL; 1656 CeedElemRestriction rstr_points = NULL; 1657 1658 CeedCallBackend(CeedOperatorAtPointsGetPoints(op, &rstr_points, &point_coords)); 1659 CeedCallBackend(CeedElemRestrictionCreateVector(rstr_points, NULL, &impl->point_coords_elem)); 1660 CeedCallBackend(CeedElemRestrictionApply(rstr_points, CEED_NOTRANSPOSE, point_coords, impl->point_coords_elem, request)); 1661 } 1662 1663 // Input Evecs and Restriction 1664 for (CeedInt i = 0; i < num_input_fields; i++) { 1665 CeedCallBackend(CeedOperatorInputRestrict_Cuda(op_input_fields[i], qf_input_fields[i], i, NULL, true, &e_data[i], impl, request)); 1666 } 1667 1668 // Clear active input Qvecs 1669 for (CeedInt i = 0; i < num_input_fields; i++) { 1670 CeedVector vec; 1671 1672 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 1673 if (vec != CEED_VECTOR_ACTIVE) continue; 1674 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 1675 } 1676 1677 // Input basis apply if needed 1678 for (CeedInt i = 0; i < num_input_fields; i++) { 1679 CeedCallBackend(CeedOperatorInputBasisAtPoints_Cuda(op_input_fields[i], qf_input_fields[i], i, num_elem, num_points, true, e_data[i], impl)); 1680 } 1681 1682 // Output pointers, as necessary 1683 for (CeedInt i = 0; i < num_output_fields; i++) { 1684 CeedEvalMode eval_mode; 1685 1686 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1687 if (eval_mode == CEED_EVAL_NONE) { 1688 // Set the output Q-Vector to use the E-Vector data directly. 1689 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields])); 1690 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields])); 1691 } 1692 } 1693 1694 // Loop over active fields 1695 for (CeedInt i = 0; i < num_input_fields; i++) { 1696 bool is_active_at_points = true; 1697 CeedInt elem_size = 1, num_comp_active = 1, e_vec_size = 0; 1698 CeedRestrictionType rstr_type; 1699 CeedVector vec; 1700 CeedElemRestriction elem_rstr; 1701 1702 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 1703 // -- Skip non-active input 1704 if (vec != CEED_VECTOR_ACTIVE) continue; 1705 1706 // -- Get active restriction type 1707 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 1708 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 1709 is_active_at_points = rstr_type == CEED_RESTRICTION_POINTS; 1710 if (!is_active_at_points) CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1711 else elem_size = max_num_points; 1712 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp_active)); 1713 1714 e_vec_size = elem_size * num_comp_active; 1715 for (CeedInt s = 0; s < e_vec_size; s++) { 1716 bool is_active_input = false; 1717 CeedEvalMode eval_mode; 1718 CeedVector vec; 1719 CeedBasis basis; 1720 1721 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 1722 // Skip non-active input 1723 is_active_input = vec == CEED_VECTOR_ACTIVE; 1724 if (!is_active_input) continue; 1725 1726 // Update unit vector 1727 if (s == 0) CeedCallBackend(CeedVectorSetValue(impl->e_vecs[i], 0.0)); 1728 else CeedCallBackend(CeedVectorSetValueStrided(impl->e_vecs[i], s - 1, e_vec_size, 0.0)); 1729 CeedCallBackend(CeedVectorSetValueStrided(impl->e_vecs[i], s, e_vec_size, 1.0)); 1730 1731 // Basis action 1732 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 1733 switch (eval_mode) { 1734 case CEED_EVAL_NONE: 1735 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i])); 1736 break; 1737 case CEED_EVAL_INTERP: 1738 case CEED_EVAL_GRAD: 1739 case CEED_EVAL_DIV: 1740 case CEED_EVAL_CURL: 1741 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 1742 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_NOTRANSPOSE, eval_mode, impl->point_coords_elem, impl->e_vecs[i], 1743 impl->q_vecs_in[i])); 1744 break; 1745 case CEED_EVAL_WEIGHT: 1746 break; // No action 1747 } 1748 1749 // Q function 1750 CeedCallBackend(CeedQFunctionApply(qf, num_elem * max_num_points, impl->q_vecs_in, impl->q_vecs_out)); 1751 1752 // Output basis apply if needed 1753 for (CeedInt j = 0; j < num_output_fields; j++) { 1754 bool is_active_output = false; 1755 CeedInt elem_size = 0; 1756 CeedRestrictionType rstr_type; 1757 CeedEvalMode eval_mode; 1758 CeedVector vec; 1759 CeedElemRestriction elem_rstr; 1760 CeedBasis basis; 1761 1762 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[j], &vec)); 1763 // ---- Skip non-active output 1764 is_active_output = vec == CEED_VECTOR_ACTIVE; 1765 if (!is_active_output) continue; 1766 1767 // ---- Check if elem size matches 1768 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr)); 1769 CeedCallBackend(CeedElemRestrictionGetType(elem_rstr, &rstr_type)); 1770 if (is_active_at_points && rstr_type != CEED_RESTRICTION_POINTS) continue; 1771 if (rstr_type == CEED_RESTRICTION_POINTS) { 1772 CeedCallBackend(CeedElemRestrictionGetMaxPointsInElement(elem_rstr, &elem_size)); 1773 } else { 1774 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 1775 } 1776 { 1777 CeedInt num_comp = 0; 1778 1779 CeedCallBackend(CeedElemRestrictionGetNumComponents(elem_rstr, &num_comp)); 1780 if (e_vec_size != num_comp * elem_size) continue; 1781 } 1782 1783 // Basis action 1784 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[j], &eval_mode)); 1785 switch (eval_mode) { 1786 case CEED_EVAL_NONE: 1787 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[j + impl->num_inputs], &e_data[j + num_input_fields])); 1788 break; 1789 case CEED_EVAL_INTERP: 1790 case CEED_EVAL_GRAD: 1791 case CEED_EVAL_DIV: 1792 case CEED_EVAL_CURL: 1793 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[j], &basis)); 1794 CeedCallBackend(CeedBasisApplyAtPoints(basis, num_elem, num_points, CEED_TRANSPOSE, eval_mode, impl->point_coords_elem, 1795 impl->q_vecs_out[j], impl->e_vecs[j + impl->num_inputs])); 1796 break; 1797 // LCOV_EXCL_START 1798 case CEED_EVAL_WEIGHT: { 1799 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 1800 // LCOV_EXCL_STOP 1801 } 1802 } 1803 1804 // Mask output e-vec 1805 CeedCallBackend(CeedVectorPointwiseMult(impl->e_vecs[j + impl->num_inputs], impl->e_vecs[i], impl->e_vecs[j + impl->num_inputs])); 1806 1807 // Restrict 1808 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[j], &elem_rstr)); 1809 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[j + impl->num_inputs], assembled, request)); 1810 1811 // Reset q_vec for 1812 if (eval_mode == CEED_EVAL_NONE) { 1813 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[j + impl->num_inputs], CEED_MEM_DEVICE, &e_data[j + num_input_fields])); 1814 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[j], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[j + num_input_fields])); 1815 } 1816 } 1817 1818 // Reset vec 1819 if (s == e_vec_size - 1 && i != num_input_fields - 1) CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 1820 } 1821 } 1822 1823 // Restore CEED_EVAL_NONE 1824 for (CeedInt i = 0; i < num_output_fields; i++) { 1825 CeedEvalMode eval_mode; 1826 1827 // Get eval_mode 1828 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1829 1830 // Restore evec 1831 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 1832 if (eval_mode == CEED_EVAL_NONE) { 1833 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields])); 1834 } 1835 } 1836 1837 // Restore input arrays 1838 for (CeedInt i = 0; i < num_input_fields; i++) { 1839 CeedCallBackend(CeedOperatorInputRestore_Cuda(op_input_fields[i], qf_input_fields[i], i, true, &e_data[i], impl)); 1840 } 1841 return CEED_ERROR_SUCCESS; 1842 } 1843 1844 //------------------------------------------------------------------------------ 1845 // Create operator 1846 //------------------------------------------------------------------------------ 1847 int CeedOperatorCreate_Cuda(CeedOperator op) { 1848 Ceed ceed; 1849 CeedOperator_Cuda *impl; 1850 1851 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1852 CeedCallBackend(CeedCalloc(1, &impl)); 1853 CeedCallBackend(CeedOperatorSetData(op, impl)); 1854 1855 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda)); 1856 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda)); 1857 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda)); 1858 CeedCallBackend( 1859 CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda)); 1860 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda)); 1861 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda)); 1862 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1863 return CEED_ERROR_SUCCESS; 1864 } 1865 1866 //------------------------------------------------------------------------------ 1867 // Create operator AtPoints 1868 //------------------------------------------------------------------------------ 1869 int CeedOperatorCreateAtPoints_Cuda(CeedOperator op) { 1870 Ceed ceed; 1871 CeedOperator_Cuda *impl; 1872 1873 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1874 CeedCallBackend(CeedCalloc(1, &impl)); 1875 CeedCallBackend(CeedOperatorSetData(op, impl)); 1876 1877 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunctionAtPoints_Cuda)); 1878 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonalAtPoints_Cuda)); 1879 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAddAtPoints_Cuda)); 1880 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1881 return CEED_ERROR_SUCCESS; 1882 } 1883 1884 //------------------------------------------------------------------------------ 1885