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