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