1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include <ceed.h> 9 #include <ceed/backend.h> 10 #include <ceed/jit-tools.h> 11 #include <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 for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) { 31 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs[i])); 32 } 33 CeedCallBackend(CeedFree(&impl->e_vecs)); 34 35 for (CeedInt i = 0; i < impl->num_inputs; i++) { 36 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i])); 37 } 38 CeedCallBackend(CeedFree(&impl->q_vecs_in)); 39 40 for (CeedInt i = 0; i < impl->num_outputs; i++) { 41 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 42 } 43 CeedCallBackend(CeedFree(&impl->q_vecs_out)); 44 45 // QFunction assembly data 46 for (CeedInt i = 0; i < impl->num_active_in; i++) { 47 CeedCallBackend(CeedVectorDestroy(&impl->qf_active_in[i])); 48 } 49 CeedCallBackend(CeedFree(&impl->qf_active_in)); 50 51 // Diag data 52 if (impl->diag) { 53 Ceed ceed; 54 55 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 56 if (impl->diag->module) { 57 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module)); 58 } 59 if (impl->diag->module_point_block) { 60 CeedCallCuda(ceed, cuModuleUnload(impl->diag->module_point_block)); 61 } 62 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_in)); 63 CeedCallCuda(ceed, cudaFree(impl->diag->d_eval_modes_out)); 64 CeedCallCuda(ceed, cudaFree(impl->diag->d_identity)); 65 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_in)); 66 CeedCallCuda(ceed, cudaFree(impl->diag->d_interp_out)); 67 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_in)); 68 CeedCallCuda(ceed, cudaFree(impl->diag->d_grad_out)); 69 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_in)); 70 CeedCallCuda(ceed, cudaFree(impl->diag->d_div_out)); 71 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_in)); 72 CeedCallCuda(ceed, cudaFree(impl->diag->d_curl_out)); 73 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->diag_rstr)); 74 CeedCallBackend(CeedElemRestrictionDestroy(&impl->diag->point_block_diag_rstr)); 75 CeedCallBackend(CeedVectorDestroy(&impl->diag->elem_diag)); 76 CeedCallBackend(CeedVectorDestroy(&impl->diag->point_block_elem_diag)); 77 } 78 CeedCallBackend(CeedFree(&impl->diag)); 79 80 if (impl->asmb) { 81 Ceed ceed; 82 83 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 84 CeedCallCuda(ceed, cuModuleUnload(impl->asmb->module)); 85 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_in)); 86 CeedCallCuda(ceed, cudaFree(impl->asmb->d_B_out)); 87 } 88 CeedCallBackend(CeedFree(&impl->asmb)); 89 90 CeedCallBackend(CeedFree(&impl)); 91 return CEED_ERROR_SUCCESS; 92 } 93 94 //------------------------------------------------------------------------------ 95 // Setup infields or outfields 96 //------------------------------------------------------------------------------ 97 static int CeedOperatorSetupFields_Cuda(CeedQFunction qf, CeedOperator op, bool is_input, CeedVector *e_vecs, CeedVector *q_vecs, CeedInt start_e, 98 CeedInt num_fields, CeedInt Q, CeedInt num_elem) { 99 Ceed ceed; 100 CeedQFunctionField *qf_fields; 101 CeedOperatorField *op_fields; 102 103 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 104 if (is_input) { 105 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 106 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 107 } else { 108 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 109 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 110 } 111 112 // Loop over fields 113 for (CeedInt i = 0; i < num_fields; i++) { 114 bool is_strided = false, skip_restriction = false; 115 CeedSize q_size; 116 CeedInt size; 117 CeedEvalMode eval_mode; 118 CeedBasis basis; 119 120 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 121 if (eval_mode != CEED_EVAL_WEIGHT) { 122 CeedElemRestriction elem_rstr; 123 124 // Check whether this field can skip the element restriction: 125 // Must be passive input, with eval_mode NONE, and have a strided restriction with CEED_STRIDES_BACKEND. 126 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &elem_rstr)); 127 128 // First, check whether the field is input or output: 129 if (is_input) { 130 CeedVector vec; 131 132 // Check for passive input 133 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 134 if (vec != CEED_VECTOR_ACTIVE) { 135 // Check eval_mode 136 if (eval_mode == CEED_EVAL_NONE) { 137 // Check for strided restriction 138 CeedCallBackend(CeedElemRestrictionIsStrided(elem_rstr, &is_strided)); 139 if (is_strided) { 140 // Check if vector is already in preferred backend ordering 141 CeedCallBackend(CeedElemRestrictionHasBackendStrides(elem_rstr, &skip_restriction)); 142 } 143 } 144 } 145 } 146 if (skip_restriction) { 147 // We do not need an E-Vector, but will use the input field vector's data directly in the operator application. 148 e_vecs[i + start_e] = NULL; 149 } else { 150 CeedCallBackend(CeedElemRestrictionCreateVector(elem_rstr, NULL, &e_vecs[i + start_e])); 151 } 152 } 153 154 switch (eval_mode) { 155 case CEED_EVAL_NONE: 156 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 157 q_size = (CeedSize)num_elem * Q * size; 158 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 159 break; 160 case CEED_EVAL_INTERP: 161 case CEED_EVAL_GRAD: 162 case CEED_EVAL_DIV: 163 case CEED_EVAL_CURL: 164 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 165 q_size = (CeedSize)num_elem * Q * size; 166 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 167 break; 168 case CEED_EVAL_WEIGHT: // Only on input fields 169 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 170 q_size = (CeedSize)num_elem * Q; 171 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 172 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i])); 173 break; 174 } 175 } 176 return CEED_ERROR_SUCCESS; 177 } 178 179 //------------------------------------------------------------------------------ 180 // CeedOperator needs to connect all the named fields (be they active or passive) to the named inputs and outputs of its CeedQFunction. 181 //------------------------------------------------------------------------------ 182 static int CeedOperatorSetup_Cuda(CeedOperator op) { 183 Ceed ceed; 184 bool is_setup_done; 185 CeedInt Q, num_elem, num_input_fields, num_output_fields; 186 CeedQFunctionField *qf_input_fields, *qf_output_fields; 187 CeedQFunction qf; 188 CeedOperatorField *op_input_fields, *op_output_fields; 189 CeedOperator_Cuda *impl; 190 191 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 192 if (is_setup_done) return CEED_ERROR_SUCCESS; 193 194 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 195 CeedCallBackend(CeedOperatorGetData(op, &impl)); 196 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 197 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 198 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 199 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 200 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 201 202 // Allocate 203 CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs)); 204 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in)); 205 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out)); 206 impl->num_inputs = num_input_fields; 207 impl->num_outputs = num_output_fields; 208 209 // Set up infield and outfield e_vecs and q_vecs 210 // Infields 211 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, true, impl->e_vecs, impl->q_vecs_in, 0, num_input_fields, Q, num_elem)); 212 // Outfields 213 CeedCallBackend(CeedOperatorSetupFields_Cuda(qf, op, false, impl->e_vecs, impl->q_vecs_out, num_input_fields, num_output_fields, Q, num_elem)); 214 215 CeedCallBackend(CeedOperatorSetSetupDone(op)); 216 return CEED_ERROR_SUCCESS; 217 } 218 219 //------------------------------------------------------------------------------ 220 // Setup Operator Inputs 221 //------------------------------------------------------------------------------ 222 static inline int CeedOperatorSetupInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 223 CeedVector in_vec, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], 224 CeedOperator_Cuda *impl, CeedRequest *request) { 225 for (CeedInt i = 0; i < num_input_fields; i++) { 226 CeedEvalMode eval_mode; 227 CeedVector vec; 228 CeedElemRestriction elem_rstr; 229 230 // Get input vector 231 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 232 if (vec == CEED_VECTOR_ACTIVE) { 233 if (skip_active) continue; 234 else vec = in_vec; 235 } 236 237 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 238 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 239 } else { 240 // Get input vector 241 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 242 // Get input element restriction 243 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 244 if (vec == CEED_VECTOR_ACTIVE) vec = in_vec; 245 // Restrict, if necessary 246 if (!impl->e_vecs[i]) { 247 // No restriction for this field; read data directly from vec. 248 CeedCallBackend(CeedVectorGetArrayRead(vec, CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i])); 249 } else { 250 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_NOTRANSPOSE, vec, impl->e_vecs[i], request)); 251 // Get evec 252 CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs[i], CEED_MEM_DEVICE, (const CeedScalar **)&e_data[i])); 253 } 254 } 255 } 256 return CEED_ERROR_SUCCESS; 257 } 258 259 //------------------------------------------------------------------------------ 260 // Input Basis Action 261 //------------------------------------------------------------------------------ 262 static inline int CeedOperatorInputBasis_Cuda(CeedInt num_elem, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 263 CeedInt num_input_fields, const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], 264 CeedOperator_Cuda *impl) { 265 for (CeedInt i = 0; i < num_input_fields; i++) { 266 CeedInt elem_size, size; 267 CeedEvalMode eval_mode; 268 CeedElemRestriction elem_rstr; 269 CeedBasis basis; 270 271 // Skip active input 272 if (skip_active) { 273 CeedVector vec; 274 275 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 276 if (vec == CEED_VECTOR_ACTIVE) continue; 277 } 278 // Get elem_size, eval_mode, size 279 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 280 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 281 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 282 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 283 // Basis action 284 switch (eval_mode) { 285 case CEED_EVAL_NONE: 286 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i])); 287 break; 288 case CEED_EVAL_INTERP: 289 case CEED_EVAL_GRAD: 290 case CEED_EVAL_DIV: 291 case CEED_EVAL_CURL: 292 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 293 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs[i], impl->q_vecs_in[i])); 294 break; 295 case CEED_EVAL_WEIGHT: 296 break; // No action 297 } 298 } 299 return CEED_ERROR_SUCCESS; 300 } 301 302 //------------------------------------------------------------------------------ 303 // Restore Input Vectors 304 //------------------------------------------------------------------------------ 305 static inline int CeedOperatorRestoreInputs_Cuda(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 306 const bool skip_active, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Cuda *impl) { 307 for (CeedInt i = 0; i < num_input_fields; i++) { 308 CeedEvalMode eval_mode; 309 CeedVector vec; 310 311 // Skip active input 312 if (skip_active) { 313 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 314 if (vec == CEED_VECTOR_ACTIVE) continue; 315 } 316 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 317 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 318 } else { 319 if (!impl->e_vecs[i]) { // This was a skip_restriction case 320 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 321 CeedCallBackend(CeedVectorRestoreArrayRead(vec, (const CeedScalar **)&e_data[i])); 322 } else { 323 CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs[i], (const CeedScalar **)&e_data[i])); 324 } 325 } 326 } 327 return CEED_ERROR_SUCCESS; 328 } 329 330 //------------------------------------------------------------------------------ 331 // Apply and add to output 332 //------------------------------------------------------------------------------ 333 static int CeedOperatorApplyAdd_Cuda(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 334 CeedInt Q, num_elem, elem_size, num_input_fields, num_output_fields, size; 335 CeedScalar *e_data[2 * CEED_FIELD_MAX] = {NULL}; 336 CeedQFunctionField *qf_input_fields, *qf_output_fields; 337 CeedQFunction qf; 338 CeedOperatorField *op_input_fields, *op_output_fields; 339 CeedOperator_Cuda *impl; 340 341 CeedCallBackend(CeedOperatorGetData(op, &impl)); 342 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 343 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 344 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 345 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 346 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 347 348 // Setup 349 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 350 351 // Input Evecs and Restriction 352 CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, in_vec, false, e_data, impl, request)); 353 354 // Input basis apply if needed 355 CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, false, e_data, impl)); 356 357 // Output pointers, as necessary 358 for (CeedInt i = 0; i < num_output_fields; i++) { 359 CeedEvalMode eval_mode; 360 361 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 362 if (eval_mode == CEED_EVAL_NONE) { 363 // Set the output Q-Vector to use the E-Vector data directly. 364 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs[i + impl->num_inputs], CEED_MEM_DEVICE, &e_data[i + num_input_fields])); 365 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_DEVICE, CEED_USE_POINTER, e_data[i + num_input_fields])); 366 } 367 } 368 369 // Q function 370 CeedCallBackend(CeedQFunctionApply(qf, num_elem * Q, impl->q_vecs_in, impl->q_vecs_out)); 371 372 // Output basis apply if needed 373 for (CeedInt i = 0; i < num_output_fields; i++) { 374 CeedEvalMode eval_mode; 375 CeedElemRestriction elem_rstr; 376 CeedBasis basis; 377 378 // Get elem_size, eval_mode, size 379 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 380 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 381 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 382 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 383 // Basis action 384 switch (eval_mode) { 385 case CEED_EVAL_NONE: 386 break; // No action 387 case CEED_EVAL_INTERP: 388 case CEED_EVAL_GRAD: 389 case CEED_EVAL_DIV: 390 case CEED_EVAL_CURL: 391 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 392 CeedCallBackend(CeedBasisApply(basis, num_elem, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs[i + impl->num_inputs])); 393 break; 394 // LCOV_EXCL_START 395 case CEED_EVAL_WEIGHT: { 396 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 397 // LCOV_EXCL_STOP 398 } 399 } 400 } 401 402 // Output restriction 403 for (CeedInt i = 0; i < num_output_fields; i++) { 404 CeedEvalMode eval_mode; 405 CeedVector vec; 406 CeedElemRestriction elem_rstr; 407 408 // Restore evec 409 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 410 if (eval_mode == CEED_EVAL_NONE) { 411 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs[i + impl->num_inputs], &e_data[i + num_input_fields])); 412 } 413 // Get output vector 414 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 415 // Restrict 416 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_output_fields[i], &elem_rstr)); 417 // Active 418 if (vec == CEED_VECTOR_ACTIVE) vec = out_vec; 419 420 CeedCallBackend(CeedElemRestrictionApply(elem_rstr, CEED_TRANSPOSE, impl->e_vecs[i + impl->num_inputs], vec, request)); 421 } 422 423 // Restore input arrays 424 CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, false, e_data, impl)); 425 return CEED_ERROR_SUCCESS; 426 } 427 428 //------------------------------------------------------------------------------ 429 // Linear QFunction Assembly Core 430 //------------------------------------------------------------------------------ 431 static inline int CeedOperatorLinearAssembleQFunctionCore_Cuda(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 432 CeedRequest *request) { 433 Ceed ceed, ceed_parent; 434 CeedInt num_active_in, num_active_out, Q, num_elem, num_input_fields, num_output_fields, size; 435 CeedScalar *assembled_array, *e_data[2 * CEED_FIELD_MAX] = {NULL}; 436 CeedVector *active_inputs; 437 CeedQFunctionField *qf_input_fields, *qf_output_fields; 438 CeedQFunction qf; 439 CeedOperatorField *op_input_fields, *op_output_fields; 440 CeedOperator_Cuda *impl; 441 442 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 443 CeedCallBackend(CeedOperatorGetFallbackParentCeed(op, &ceed_parent)); 444 CeedCallBackend(CeedOperatorGetData(op, &impl)); 445 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 446 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 447 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 448 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 449 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 450 active_inputs = impl->qf_active_in; 451 num_active_in = impl->num_active_in, num_active_out = impl->num_active_out; 452 453 // Setup 454 CeedCallBackend(CeedOperatorSetup_Cuda(op)); 455 456 // Input Evecs and Restriction 457 CeedCallBackend(CeedOperatorSetupInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, NULL, true, e_data, impl, request)); 458 459 // Count number of active input fields 460 if (!num_active_in) { 461 for (CeedInt i = 0; i < num_input_fields; i++) { 462 CeedScalar *q_vec_array; 463 CeedVector vec; 464 465 // Get input vector 466 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 467 // Check if active input 468 if (vec == CEED_VECTOR_ACTIVE) { 469 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 470 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 471 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_DEVICE, &q_vec_array)); 472 CeedCallBackend(CeedRealloc(num_active_in + size, &active_inputs)); 473 for (CeedInt field = 0; field < size; field++) { 474 CeedSize q_size = (CeedSize)Q * num_elem; 475 476 CeedCallBackend(CeedVectorCreate(ceed, q_size, &active_inputs[num_active_in + field])); 477 CeedCallBackend( 478 CeedVectorSetArray(active_inputs[num_active_in + field], CEED_MEM_DEVICE, CEED_USE_POINTER, &q_vec_array[field * Q * num_elem])); 479 } 480 num_active_in += size; 481 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &q_vec_array)); 482 } 483 } 484 impl->num_active_in = num_active_in; 485 impl->qf_active_in = active_inputs; 486 } 487 488 // Count number of active output fields 489 if (!num_active_out) { 490 for (CeedInt i = 0; i < num_output_fields; i++) { 491 CeedVector vec; 492 493 // Get output vector 494 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 495 // Check if active output 496 if (vec == CEED_VECTOR_ACTIVE) { 497 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &size)); 498 num_active_out += size; 499 } 500 } 501 impl->num_active_out = num_active_out; 502 } 503 504 // Check sizes 505 CeedCheck(num_active_in > 0 && num_active_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 506 507 // Build objects if needed 508 if (build_objects) { 509 CeedSize l_size = (CeedSize)num_elem * Q * num_active_in * num_active_out; 510 CeedInt strides[3] = {1, num_elem * Q, Q}; /* *NOPAD* */ 511 512 // Create output restriction 513 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed_parent, num_elem, Q, num_active_in * num_active_out, 514 num_active_in * num_active_out * num_elem * Q, strides, rstr)); 515 // Create assembled vector 516 CeedCallBackend(CeedVectorCreate(ceed_parent, l_size, assembled)); 517 } 518 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0)); 519 CeedCallBackend(CeedVectorGetArray(*assembled, CEED_MEM_DEVICE, &assembled_array)); 520 521 // Input basis apply 522 CeedCallBackend(CeedOperatorInputBasis_Cuda(num_elem, qf_input_fields, op_input_fields, num_input_fields, true, e_data, impl)); 523 524 // Assemble QFunction 525 for (CeedInt in = 0; in < num_active_in; in++) { 526 // Set Inputs 527 CeedCallBackend(CeedVectorSetValue(active_inputs[in], 1.0)); 528 if (num_active_in > 1) { 529 CeedCallBackend(CeedVectorSetValue(active_inputs[(in + num_active_in - 1) % num_active_in], 0.0)); 530 } 531 // Set Outputs 532 for (CeedInt out = 0; out < num_output_fields; out++) { 533 CeedVector vec; 534 535 // Get output vector 536 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 537 // Check if active output 538 if (vec == CEED_VECTOR_ACTIVE) { 539 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, CEED_USE_POINTER, assembled_array)); 540 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &size)); 541 assembled_array += size * Q * num_elem; // Advance the pointer by the size of the output 542 } 543 } 544 // Apply QFunction 545 CeedCallBackend(CeedQFunctionApply(qf, Q * num_elem, impl->q_vecs_in, impl->q_vecs_out)); 546 } 547 548 // Un-set output q_vecs to prevent accidental overwrite of Assembled 549 for (CeedInt out = 0; out < num_output_fields; out++) { 550 CeedVector vec; 551 552 // Get output vector 553 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 554 // Check if active output 555 if (vec == CEED_VECTOR_ACTIVE) { 556 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_DEVICE, NULL)); 557 } 558 } 559 560 // Restore input arrays 561 CeedCallBackend(CeedOperatorRestoreInputs_Cuda(num_input_fields, qf_input_fields, op_input_fields, true, e_data, impl)); 562 563 // Restore output 564 CeedCallBackend(CeedVectorRestoreArray(*assembled, &assembled_array)); 565 return CEED_ERROR_SUCCESS; 566 } 567 568 //------------------------------------------------------------------------------ 569 // Assemble Linear QFunction 570 //------------------------------------------------------------------------------ 571 static int CeedOperatorLinearAssembleQFunction_Cuda(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 572 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, true, assembled, rstr, request); 573 } 574 575 //------------------------------------------------------------------------------ 576 // Update Assembled Linear QFunction 577 //------------------------------------------------------------------------------ 578 static int CeedOperatorLinearAssembleQFunctionUpdate_Cuda(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 579 return CeedOperatorLinearAssembleQFunctionCore_Cuda(op, false, &assembled, &rstr, request); 580 } 581 582 //------------------------------------------------------------------------------ 583 // Assemble Diagonal Setup 584 //------------------------------------------------------------------------------ 585 static inline int CeedOperatorAssembleDiagonalSetup_Cuda(CeedOperator op) { 586 Ceed ceed; 587 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 588 CeedInt q_comp, num_nodes, num_qpts; 589 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 590 CeedBasis basis_in = NULL, basis_out = NULL; 591 CeedQFunctionField *qf_fields; 592 CeedQFunction qf; 593 CeedOperatorField *op_fields; 594 CeedOperator_Cuda *impl; 595 596 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 597 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 598 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 599 600 // Determine active input basis 601 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 602 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 603 for (CeedInt i = 0; i < num_input_fields; i++) { 604 CeedVector vec; 605 606 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 607 if (vec == CEED_VECTOR_ACTIVE) { 608 CeedBasis basis; 609 CeedEvalMode eval_mode; 610 611 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 612 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, 613 "Backend does not implement operator diagonal assembly with multiple active bases"); 614 basis_in = basis; 615 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 616 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 617 if (eval_mode != CEED_EVAL_WEIGHT) { 618 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 619 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 620 for (CeedInt d = 0; d < q_comp; d++) eval_modes_in[num_eval_modes_in + d] = eval_mode; 621 num_eval_modes_in += q_comp; 622 } 623 } 624 } 625 626 // Determine active output basis 627 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 628 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 629 for (CeedInt i = 0; i < num_output_fields; i++) { 630 CeedVector vec; 631 632 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 633 if (vec == CEED_VECTOR_ACTIVE) { 634 CeedBasis basis; 635 CeedEvalMode eval_mode; 636 637 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 638 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 639 "Backend does not implement operator diagonal assembly with multiple active bases"); 640 basis_out = basis; 641 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 642 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 643 if (eval_mode != CEED_EVAL_WEIGHT) { 644 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF assembly 645 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 646 for (CeedInt d = 0; d < q_comp; d++) eval_modes_out[num_eval_modes_out + d] = eval_mode; 647 num_eval_modes_out += q_comp; 648 } 649 } 650 } 651 652 // Operator data struct 653 CeedCallBackend(CeedOperatorGetData(op, &impl)); 654 CeedCallBackend(CeedCalloc(1, &impl->diag)); 655 CeedOperatorDiag_Cuda *diag = impl->diag; 656 657 // Basis matrices 658 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 659 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 660 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 661 const CeedInt interp_bytes = num_nodes * num_qpts * sizeof(CeedScalar); 662 const CeedInt eval_modes_bytes = sizeof(CeedEvalMode); 663 bool has_eval_none = false; 664 665 // CEED_EVAL_NONE 666 for (CeedInt i = 0; i < num_eval_modes_in; i++) has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 667 for (CeedInt i = 0; i < num_eval_modes_out; i++) has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 668 if (has_eval_none) { 669 CeedScalar *identity = NULL; 670 671 CeedCallBackend(CeedCalloc(num_nodes * num_qpts, &identity)); 672 for (CeedInt i = 0; i < (num_nodes < num_qpts ? num_nodes : num_qpts); i++) identity[i * num_nodes + i] = 1.0; 673 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_identity, interp_bytes)); 674 CeedCallCuda(ceed, cudaMemcpy(diag->d_identity, identity, interp_bytes, cudaMemcpyHostToDevice)); 675 CeedCallBackend(CeedFree(&identity)); 676 } 677 678 // CEED_EVAL_INTERP, CEED_EVAL_GRAD, CEED_EVAL_DIV, and CEED_EVAL_CURL 679 for (CeedInt in = 0; in < 2; in++) { 680 CeedFESpace fespace; 681 CeedBasis basis = in ? basis_in : basis_out; 682 683 CeedCallBackend(CeedBasisGetFESpace(basis, &fespace)); 684 switch (fespace) { 685 case CEED_FE_SPACE_H1: { 686 CeedInt q_comp_interp, q_comp_grad; 687 const CeedScalar *interp, *grad; 688 CeedScalar *d_interp, *d_grad; 689 690 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 691 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_GRAD, &q_comp_grad)); 692 693 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 694 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 695 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 696 CeedCallBackend(CeedBasisGetGrad(basis, &grad)); 697 CeedCallCuda(ceed, cudaMalloc((void **)&d_grad, interp_bytes * q_comp_grad)); 698 CeedCallCuda(ceed, cudaMemcpy(d_grad, grad, interp_bytes * q_comp_grad, cudaMemcpyHostToDevice)); 699 if (in) { 700 diag->d_interp_in = d_interp; 701 diag->d_grad_in = d_grad; 702 } else { 703 diag->d_interp_out = d_interp; 704 diag->d_grad_out = d_grad; 705 } 706 } break; 707 case CEED_FE_SPACE_HDIV: { 708 CeedInt q_comp_interp, q_comp_div; 709 const CeedScalar *interp, *div; 710 CeedScalar *d_interp, *d_div; 711 712 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 713 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_DIV, &q_comp_div)); 714 715 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 716 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 717 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 718 CeedCallBackend(CeedBasisGetDiv(basis, &div)); 719 CeedCallCuda(ceed, cudaMalloc((void **)&d_div, interp_bytes * q_comp_div)); 720 CeedCallCuda(ceed, cudaMemcpy(d_div, div, interp_bytes * q_comp_div, cudaMemcpyHostToDevice)); 721 if (in) { 722 diag->d_interp_in = d_interp; 723 diag->d_div_in = d_div; 724 } else { 725 diag->d_interp_out = d_interp; 726 diag->d_div_out = d_div; 727 } 728 } break; 729 case CEED_FE_SPACE_HCURL: { 730 CeedInt q_comp_interp, q_comp_curl; 731 const CeedScalar *interp, *curl; 732 CeedScalar *d_interp, *d_curl; 733 734 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_INTERP, &q_comp_interp)); 735 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis, CEED_EVAL_CURL, &q_comp_curl)); 736 737 CeedCallBackend(CeedBasisGetInterp(basis, &interp)); 738 CeedCallCuda(ceed, cudaMalloc((void **)&d_interp, interp_bytes * q_comp_interp)); 739 CeedCallCuda(ceed, cudaMemcpy(d_interp, interp, interp_bytes * q_comp_interp, cudaMemcpyHostToDevice)); 740 CeedCallBackend(CeedBasisGetCurl(basis, &curl)); 741 CeedCallCuda(ceed, cudaMalloc((void **)&d_curl, interp_bytes * q_comp_curl)); 742 CeedCallCuda(ceed, cudaMemcpy(d_curl, curl, interp_bytes * q_comp_curl, cudaMemcpyHostToDevice)); 743 if (in) { 744 diag->d_interp_in = d_interp; 745 diag->d_curl_in = d_curl; 746 } else { 747 diag->d_interp_out = d_interp; 748 diag->d_curl_out = d_curl; 749 } 750 } break; 751 } 752 } 753 754 // Arrays of eval_modes 755 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_in, num_eval_modes_in * eval_modes_bytes)); 756 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_in, eval_modes_in, num_eval_modes_in * eval_modes_bytes, cudaMemcpyHostToDevice)); 757 CeedCallCuda(ceed, cudaMalloc((void **)&diag->d_eval_modes_out, num_eval_modes_out * eval_modes_bytes)); 758 CeedCallCuda(ceed, cudaMemcpy(diag->d_eval_modes_out, eval_modes_out, num_eval_modes_out * eval_modes_bytes, cudaMemcpyHostToDevice)); 759 CeedCallBackend(CeedFree(&eval_modes_in)); 760 CeedCallBackend(CeedFree(&eval_modes_out)); 761 return CEED_ERROR_SUCCESS; 762 } 763 764 //------------------------------------------------------------------------------ 765 // Assemble Diagonal Setup (Compilation) 766 //------------------------------------------------------------------------------ 767 static inline int CeedOperatorAssembleDiagonalSetupCompile_Cuda(CeedOperator op, CeedInt use_ceedsize_idx, const bool is_point_block) { 768 Ceed ceed; 769 char *diagonal_kernel_source; 770 const char *diagonal_kernel_path; 771 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 772 CeedInt num_comp, q_comp, num_nodes, num_qpts; 773 CeedBasis basis_in = NULL, basis_out = NULL; 774 CeedQFunctionField *qf_fields; 775 CeedQFunction qf; 776 CeedOperatorField *op_fields; 777 CeedOperator_Cuda *impl; 778 779 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 780 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 781 CeedCallBackend(CeedQFunctionGetNumArgs(qf, &num_input_fields, &num_output_fields)); 782 783 // Determine active input basis 784 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 785 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 786 for (CeedInt i = 0; i < num_input_fields; i++) { 787 CeedVector vec; 788 789 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 790 if (vec == CEED_VECTOR_ACTIVE) { 791 CeedEvalMode eval_mode; 792 793 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_in)); 794 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 795 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 796 if (eval_mode != CEED_EVAL_WEIGHT) { 797 num_eval_modes_in += q_comp; 798 } 799 } 800 } 801 802 // Determine active output basis 803 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 804 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 805 for (CeedInt i = 0; i < num_output_fields; i++) { 806 CeedVector vec; 807 808 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec)); 809 if (vec == CEED_VECTOR_ACTIVE) { 810 CeedEvalMode eval_mode; 811 812 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis_out)); 813 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 814 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 815 if (eval_mode != CEED_EVAL_WEIGHT) { 816 num_eval_modes_out += q_comp; 817 } 818 } 819 } 820 821 // Operator data struct 822 CeedCallBackend(CeedOperatorGetData(op, &impl)); 823 CeedOperatorDiag_Cuda *diag = impl->diag; 824 825 // Assemble kernel 826 CUmodule *module = is_point_block ? &diag->module_point_block : &diag->module; 827 CeedInt elems_per_block = 1; 828 CeedCallBackend(CeedBasisGetNumNodes(basis_in, &num_nodes)); 829 CeedCallBackend(CeedBasisGetNumComponents(basis_in, &num_comp)); 830 if (basis_in == CEED_BASIS_NONE) num_qpts = num_nodes; 831 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts)); 832 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble-diagonal.h", &diagonal_kernel_path)); 833 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Kernel Source -----\n"); 834 CeedCallBackend(CeedLoadSourceToBuffer(ceed, diagonal_kernel_path, &diagonal_kernel_source)); 835 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Diagonal Assembly Source Complete! -----\n"); 836 CeedCallCuda(ceed, CeedCompile_Cuda(ceed, diagonal_kernel_source, module, 8, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 837 num_eval_modes_out, "NUM_COMP", num_comp, "NUM_NODES", num_nodes, "NUM_QPTS", num_qpts, "USE_CEEDSIZE", 838 use_ceedsize_idx, "USE_POINT_BLOCK", is_point_block ? 1 : 0, "BLOCK_SIZE", num_nodes * elems_per_block)); 839 CeedCallCuda(ceed, CeedGetKernel_Cuda(ceed, *module, "LinearDiagonal", is_point_block ? &diag->LinearPointBlock : &diag->LinearDiagonal)); 840 CeedCallBackend(CeedFree(&diagonal_kernel_path)); 841 CeedCallBackend(CeedFree(&diagonal_kernel_source)); 842 return CEED_ERROR_SUCCESS; 843 } 844 845 //------------------------------------------------------------------------------ 846 // Assemble Diagonal Core 847 //------------------------------------------------------------------------------ 848 static inline int CeedOperatorAssembleDiagonalCore_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request, const bool is_point_block) { 849 Ceed ceed; 850 CeedInt num_elem, num_nodes; 851 CeedScalar *elem_diag_array; 852 const CeedScalar *assembled_qf_array; 853 CeedVector assembled_qf = NULL, elem_diag; 854 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out, diag_rstr; 855 CeedOperator_Cuda *impl; 856 857 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 858 CeedCallBackend(CeedOperatorGetData(op, &impl)); 859 860 // Assemble QFunction 861 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, request)); 862 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 863 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 864 865 // Setup 866 if (!impl->diag) CeedCallBackend(CeedOperatorAssembleDiagonalSetup_Cuda(op)); 867 CeedOperatorDiag_Cuda *diag = impl->diag; 868 869 assert(diag != NULL); 870 871 // Assemble kernel if needed 872 if ((!is_point_block && !diag->LinearDiagonal) || (is_point_block && !diag->LinearPointBlock)) { 873 CeedSize assembled_length, assembled_qf_length; 874 CeedInt use_ceedsize_idx = 0; 875 CeedCallBackend(CeedVectorGetLength(assembled, &assembled_length)); 876 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 877 if ((assembled_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 878 879 CeedCallBackend(CeedOperatorAssembleDiagonalSetupCompile_Cuda(op, use_ceedsize_idx, is_point_block)); 880 } 881 882 // Restriction and diagonal vector 883 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 884 CeedCheck(rstr_in == rstr_out, ceed, CEED_ERROR_BACKEND, 885 "Cannot assemble operator diagonal with different input and output active element restrictions"); 886 if (!is_point_block && !diag->diag_rstr) { 887 CeedCallBackend(CeedElemRestrictionCreateUnsignedCopy(rstr_out, &diag->diag_rstr)); 888 CeedCallBackend(CeedElemRestrictionCreateVector(diag->diag_rstr, NULL, &diag->elem_diag)); 889 } else if (is_point_block && !diag->point_block_diag_rstr) { 890 CeedCallBackend(CeedOperatorCreateActivePointBlockRestriction(rstr_out, &diag->point_block_diag_rstr)); 891 CeedCallBackend(CeedElemRestrictionCreateVector(diag->point_block_diag_rstr, NULL, &diag->point_block_elem_diag)); 892 } 893 diag_rstr = is_point_block ? diag->point_block_diag_rstr : diag->diag_rstr; 894 elem_diag = is_point_block ? diag->point_block_elem_diag : diag->elem_diag; 895 CeedCallBackend(CeedVectorSetValue(elem_diag, 0.0)); 896 897 // Only assemble diagonal if the basis has nodes, otherwise inputs are null pointers 898 CeedCallBackend(CeedElemRestrictionGetElementSize(diag_rstr, &num_nodes)); 899 if (num_nodes > 0) { 900 // Assemble element operator diagonals 901 CeedCallBackend(CeedElemRestrictionGetNumElements(diag_rstr, &num_elem)); 902 CeedCallBackend(CeedVectorGetArray(elem_diag, CEED_MEM_DEVICE, &elem_diag_array)); 903 904 // Compute the diagonal of B^T D B 905 CeedInt elems_per_block = 1; 906 CeedInt grid = CeedDivUpInt(num_elem, elems_per_block); 907 void *args[] = {(void *)&num_elem, &diag->d_identity, &diag->d_interp_in, &diag->d_grad_in, &diag->d_div_in, 908 &diag->d_curl_in, &diag->d_interp_out, &diag->d_grad_out, &diag->d_div_out, &diag->d_curl_out, 909 &diag->d_eval_modes_in, &diag->d_eval_modes_out, &assembled_qf_array, &elem_diag_array}; 910 911 if (is_point_block) { 912 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearPointBlock, grid, num_nodes, 1, elems_per_block, args)); 913 } else { 914 CeedCallBackend(CeedRunKernelDim_Cuda(ceed, diag->LinearDiagonal, grid, num_nodes, 1, elems_per_block, args)); 915 } 916 917 // Restore arrays 918 CeedCallBackend(CeedVectorRestoreArray(elem_diag, &elem_diag_array)); 919 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 920 } 921 922 // Assemble local operator diagonal 923 CeedCallBackend(CeedElemRestrictionApply(diag_rstr, CEED_TRANSPOSE, elem_diag, assembled, request)); 924 925 // Cleanup 926 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 927 return CEED_ERROR_SUCCESS; 928 } 929 930 //------------------------------------------------------------------------------ 931 // Assemble Linear Diagonal 932 //------------------------------------------------------------------------------ 933 static int CeedOperatorLinearAssembleAddDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 934 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, false)); 935 return CEED_ERROR_SUCCESS; 936 } 937 938 //------------------------------------------------------------------------------ 939 // Assemble Linear Point Block Diagonal 940 //------------------------------------------------------------------------------ 941 static int CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda(CeedOperator op, CeedVector assembled, CeedRequest *request) { 942 CeedCallBackend(CeedOperatorAssembleDiagonalCore_Cuda(op, assembled, request, true)); 943 return CEED_ERROR_SUCCESS; 944 } 945 946 //------------------------------------------------------------------------------ 947 // Single Operator Assembly Setup 948 //------------------------------------------------------------------------------ 949 static int CeedSingleOperatorAssembleSetup_Cuda(CeedOperator op, CeedInt use_ceedsize_idx) { 950 Ceed ceed; 951 Ceed_Cuda *cuda_data; 952 char *assembly_kernel_source; 953 const char *assembly_kernel_path; 954 CeedInt num_input_fields, num_output_fields, num_eval_modes_in = 0, num_eval_modes_out = 0; 955 CeedInt elem_size_in, num_qpts_in = 0, num_comp_in, elem_size_out, num_qpts_out, num_comp_out, q_comp; 956 CeedEvalMode *eval_modes_in = NULL, *eval_modes_out = NULL; 957 CeedElemRestriction rstr_in = NULL, rstr_out = NULL; 958 CeedBasis basis_in = NULL, basis_out = NULL; 959 CeedQFunctionField *qf_fields; 960 CeedQFunction qf; 961 CeedOperatorField *input_fields, *output_fields; 962 CeedOperator_Cuda *impl; 963 964 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 965 CeedCallBackend(CeedOperatorGetData(op, &impl)); 966 967 // Get intput and output fields 968 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &input_fields, &num_output_fields, &output_fields)); 969 970 // Determine active input basis eval mode 971 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 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(input_fields[i], &vec)); 977 if (vec == CEED_VECTOR_ACTIVE) { 978 CeedBasis basis; 979 CeedEvalMode eval_mode; 980 981 CeedCallBackend(CeedOperatorFieldGetBasis(input_fields[i], &basis)); 982 CeedCheck(!basis_in || basis_in == basis, ceed, CEED_ERROR_BACKEND, "Backend does not implement operator assembly with multiple active bases"); 983 basis_in = basis; 984 CeedCallBackend(CeedOperatorFieldGetElemRestriction(input_fields[i], &rstr_in)); 985 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 986 if (basis_in == CEED_BASIS_NONE) num_qpts_in = elem_size_in; 987 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_in, &num_qpts_in)); 988 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 989 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_mode, &q_comp)); 990 if (eval_mode != CEED_EVAL_WEIGHT) { 991 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 992 CeedCallBackend(CeedRealloc(num_eval_modes_in + q_comp, &eval_modes_in)); 993 for (CeedInt d = 0; d < q_comp; d++) { 994 eval_modes_in[num_eval_modes_in + d] = eval_mode; 995 } 996 num_eval_modes_in += q_comp; 997 } 998 } 999 } 1000 1001 // Determine active output basis; basis_out and rstr_out only used if same as input, TODO 1002 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 1003 for (CeedInt i = 0; i < num_output_fields; i++) { 1004 CeedVector vec; 1005 1006 CeedCallBackend(CeedOperatorFieldGetVector(output_fields[i], &vec)); 1007 if (vec == CEED_VECTOR_ACTIVE) { 1008 CeedBasis basis; 1009 CeedEvalMode eval_mode; 1010 1011 CeedCallBackend(CeedOperatorFieldGetBasis(output_fields[i], &basis)); 1012 CeedCheck(!basis_out || basis_out == basis, ceed, CEED_ERROR_BACKEND, 1013 "Backend does not implement operator assembly with multiple active bases"); 1014 basis_out = basis; 1015 CeedCallBackend(CeedOperatorFieldGetElemRestriction(output_fields[i], &rstr_out)); 1016 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1017 if (basis_out == CEED_BASIS_NONE) num_qpts_out = elem_size_out; 1018 else CeedCallBackend(CeedBasisGetNumQuadraturePoints(basis_out, &num_qpts_out)); 1019 CeedCheck(num_qpts_in == num_qpts_out, ceed, CEED_ERROR_UNSUPPORTED, 1020 "Active input and output bases must have the same number of quadrature points"); 1021 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 1022 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_mode, &q_comp)); 1023 if (eval_mode != CEED_EVAL_WEIGHT) { 1024 // q_comp = 1 if CEED_EVAL_NONE, CEED_EVAL_WEIGHT caught by QF Assembly 1025 CeedCallBackend(CeedRealloc(num_eval_modes_out + q_comp, &eval_modes_out)); 1026 for (CeedInt d = 0; d < q_comp; d++) { 1027 eval_modes_out[num_eval_modes_out + d] = eval_mode; 1028 } 1029 num_eval_modes_out += q_comp; 1030 } 1031 } 1032 } 1033 CeedCheck(num_eval_modes_in > 0 && num_eval_modes_out > 0, ceed, CEED_ERROR_UNSUPPORTED, "Cannot assemble operator without inputs/outputs"); 1034 1035 CeedCallBackend(CeedCalloc(1, &impl->asmb)); 1036 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1037 asmb->elems_per_block = 1; 1038 asmb->block_size_x = elem_size_in; 1039 asmb->block_size_y = elem_size_out; 1040 1041 CeedCallBackend(CeedGetData(ceed, &cuda_data)); 1042 bool fallback = asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block > cuda_data->device_prop.maxThreadsPerBlock; 1043 1044 if (fallback) { 1045 // Use fallback kernel with 1D threadblock 1046 asmb->block_size_y = 1; 1047 } 1048 1049 // Compile kernels 1050 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_in, &num_comp_in)); 1051 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr_out, &num_comp_out)); 1052 CeedCallBackend(CeedGetJitAbsolutePath(ceed, "ceed/jit-source/cuda/cuda-ref-operator-assemble.h", &assembly_kernel_path)); 1053 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Kernel Source -----\n"); 1054 CeedCallBackend(CeedLoadSourceToBuffer(ceed, assembly_kernel_path, &assembly_kernel_source)); 1055 CeedDebug256(ceed, CEED_DEBUG_COLOR_SUCCESS, "----- Loading Assembly Source Complete! -----\n"); 1056 CeedCallBackend(CeedCompile_Cuda(ceed, assembly_kernel_source, &asmb->module, 10, "NUM_EVAL_MODES_IN", num_eval_modes_in, "NUM_EVAL_MODES_OUT", 1057 num_eval_modes_out, "NUM_COMP_IN", num_comp_in, "NUM_COMP_OUT", num_comp_out, "NUM_NODES_IN", elem_size_in, 1058 "NUM_NODES_OUT", elem_size_out, "NUM_QPTS", num_qpts_in, "BLOCK_SIZE", 1059 asmb->block_size_x * asmb->block_size_y * asmb->elems_per_block, "BLOCK_SIZE_Y", asmb->block_size_y, 1060 "USE_CEEDSIZE", use_ceedsize_idx)); 1061 CeedCallBackend(CeedGetKernel_Cuda(ceed, asmb->module, "LinearAssemble", &asmb->LinearAssemble)); 1062 CeedCallBackend(CeedFree(&assembly_kernel_path)); 1063 CeedCallBackend(CeedFree(&assembly_kernel_source)); 1064 1065 // Load into B_in, in order that they will be used in eval_modes_in 1066 { 1067 const CeedInt in_bytes = elem_size_in * num_qpts_in * num_eval_modes_in * sizeof(CeedScalar); 1068 CeedInt d_in = 0; 1069 CeedEvalMode eval_modes_in_prev = CEED_EVAL_NONE; 1070 bool has_eval_none = false; 1071 CeedScalar *identity = NULL; 1072 1073 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1074 has_eval_none = has_eval_none || (eval_modes_in[i] == CEED_EVAL_NONE); 1075 } 1076 if (has_eval_none) { 1077 CeedCallBackend(CeedCalloc(elem_size_in * num_qpts_in, &identity)); 1078 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; 1079 } 1080 1081 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_in, in_bytes)); 1082 for (CeedInt i = 0; i < num_eval_modes_in; i++) { 1083 const CeedScalar *h_B_in; 1084 1085 CeedCallBackend(CeedOperatorGetBasisPointer(basis_in, eval_modes_in[i], identity, &h_B_in)); 1086 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_in, eval_modes_in[i], &q_comp)); 1087 if (q_comp > 1) { 1088 if (i == 0 || eval_modes_in[i] != eval_modes_in_prev) d_in = 0; 1089 else h_B_in = &h_B_in[(++d_in) * elem_size_in * num_qpts_in]; 1090 } 1091 eval_modes_in_prev = eval_modes_in[i]; 1092 1093 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), 1094 cudaMemcpyHostToDevice)); 1095 } 1096 1097 if (identity) { 1098 CeedCallBackend(CeedFree(&identity)); 1099 } 1100 } 1101 1102 // Load into B_out, in order that they will be used in eval_modes_out 1103 { 1104 const CeedInt out_bytes = elem_size_out * num_qpts_out * num_eval_modes_out * sizeof(CeedScalar); 1105 CeedInt d_out = 0; 1106 CeedEvalMode eval_modes_out_prev = CEED_EVAL_NONE; 1107 bool has_eval_none = false; 1108 CeedScalar *identity = NULL; 1109 1110 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1111 has_eval_none = has_eval_none || (eval_modes_out[i] == CEED_EVAL_NONE); 1112 } 1113 if (has_eval_none) { 1114 CeedCallBackend(CeedCalloc(elem_size_out * num_qpts_out, &identity)); 1115 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; 1116 } 1117 1118 CeedCallCuda(ceed, cudaMalloc((void **)&asmb->d_B_out, out_bytes)); 1119 for (CeedInt i = 0; i < num_eval_modes_out; i++) { 1120 const CeedScalar *h_B_out; 1121 1122 CeedCallBackend(CeedOperatorGetBasisPointer(basis_out, eval_modes_out[i], identity, &h_B_out)); 1123 CeedCallBackend(CeedBasisGetNumQuadratureComponents(basis_out, eval_modes_out[i], &q_comp)); 1124 if (q_comp > 1) { 1125 if (i == 0 || eval_modes_out[i] != eval_modes_out_prev) d_out = 0; 1126 else h_B_out = &h_B_out[(++d_out) * elem_size_out * num_qpts_out]; 1127 } 1128 eval_modes_out_prev = eval_modes_out[i]; 1129 1130 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), 1131 cudaMemcpyHostToDevice)); 1132 } 1133 1134 if (identity) { 1135 CeedCallBackend(CeedFree(&identity)); 1136 } 1137 } 1138 return CEED_ERROR_SUCCESS; 1139 } 1140 1141 //------------------------------------------------------------------------------ 1142 // Assemble matrix data for COO matrix of assembled operator. 1143 // The sparsity pattern is set by CeedOperatorLinearAssembleSymbolic. 1144 // 1145 // Note that this (and other assembly routines) currently assume only one active input restriction/basis per operator (could have multiple basis eval 1146 // modes). 1147 // TODO: allow multiple active input restrictions/basis objects 1148 //------------------------------------------------------------------------------ 1149 static int CeedSingleOperatorAssemble_Cuda(CeedOperator op, CeedInt offset, CeedVector values) { 1150 Ceed ceed; 1151 CeedSize values_length = 0, assembled_qf_length = 0; 1152 CeedInt use_ceedsize_idx = 0, num_elem_in, num_elem_out, elem_size_in, elem_size_out; 1153 CeedScalar *values_array; 1154 const CeedScalar *assembled_qf_array; 1155 CeedVector assembled_qf = NULL; 1156 CeedElemRestriction assembled_rstr = NULL, rstr_in, rstr_out; 1157 CeedRestrictionType rstr_type_in, rstr_type_out; 1158 const bool *orients_in = NULL, *orients_out = NULL; 1159 const CeedInt8 *curl_orients_in = NULL, *curl_orients_out = NULL; 1160 CeedOperator_Cuda *impl; 1161 1162 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1163 CeedCallBackend(CeedOperatorGetData(op, &impl)); 1164 1165 // Assemble QFunction 1166 CeedCallBackend(CeedOperatorLinearAssembleQFunctionBuildOrUpdate(op, &assembled_qf, &assembled_rstr, CEED_REQUEST_IMMEDIATE)); 1167 CeedCallBackend(CeedElemRestrictionDestroy(&assembled_rstr)); 1168 CeedCallBackend(CeedVectorGetArrayRead(assembled_qf, CEED_MEM_DEVICE, &assembled_qf_array)); 1169 1170 CeedCallBackend(CeedVectorGetLength(values, &values_length)); 1171 CeedCallBackend(CeedVectorGetLength(assembled_qf, &assembled_qf_length)); 1172 if ((values_length > INT_MAX) || (assembled_qf_length > INT_MAX)) use_ceedsize_idx = 1; 1173 1174 // Setup 1175 if (!impl->asmb) CeedCallBackend(CeedSingleOperatorAssembleSetup_Cuda(op, use_ceedsize_idx)); 1176 CeedOperatorAssemble_Cuda *asmb = impl->asmb; 1177 1178 assert(asmb != NULL); 1179 1180 // Assemble element operator 1181 CeedCallBackend(CeedVectorGetArray(values, CEED_MEM_DEVICE, &values_array)); 1182 values_array += offset; 1183 1184 CeedCallBackend(CeedOperatorGetActiveElemRestrictions(op, &rstr_in, &rstr_out)); 1185 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_in, &num_elem_in)); 1186 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_in, &elem_size_in)); 1187 1188 CeedCallBackend(CeedElemRestrictionGetType(rstr_in, &rstr_type_in)); 1189 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1190 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_in, CEED_MEM_DEVICE, &orients_in)); 1191 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1192 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_in, CEED_MEM_DEVICE, &curl_orients_in)); 1193 } 1194 1195 if (rstr_in != rstr_out) { 1196 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr_out, &num_elem_out)); 1197 CeedCheck(num_elem_in == num_elem_out, ceed, CEED_ERROR_UNSUPPORTED, 1198 "Active input and output operator restrictions must have the same number of elements"); 1199 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr_out, &elem_size_out)); 1200 1201 CeedCallBackend(CeedElemRestrictionGetType(rstr_out, &rstr_type_out)); 1202 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1203 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr_out, CEED_MEM_DEVICE, &orients_out)); 1204 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1205 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr_out, CEED_MEM_DEVICE, &curl_orients_out)); 1206 } 1207 } else { 1208 elem_size_out = elem_size_in; 1209 orients_out = orients_in; 1210 curl_orients_out = curl_orients_in; 1211 } 1212 1213 // Compute B^T D B 1214 CeedInt shared_mem = 1215 ((curl_orients_in || curl_orients_out ? elem_size_in * elem_size_out : 0) + (curl_orients_in ? elem_size_in * asmb->block_size_y : 0)) * 1216 sizeof(CeedScalar); 1217 CeedInt grid = CeedDivUpInt(num_elem_in, asmb->elems_per_block); 1218 void *args[] = {(void *)&num_elem_in, &asmb->d_B_in, &asmb->d_B_out, &orients_in, &curl_orients_in, 1219 &orients_out, &curl_orients_out, &assembled_qf_array, &values_array}; 1220 1221 CeedCallBackend( 1222 CeedRunKernelDimShared_Cuda(ceed, asmb->LinearAssemble, grid, asmb->block_size_x, asmb->block_size_y, asmb->elems_per_block, shared_mem, args)); 1223 1224 // Restore arrays 1225 CeedCallBackend(CeedVectorRestoreArray(values, &values_array)); 1226 CeedCallBackend(CeedVectorRestoreArrayRead(assembled_qf, &assembled_qf_array)); 1227 1228 // Cleanup 1229 CeedCallBackend(CeedVectorDestroy(&assembled_qf)); 1230 if (rstr_type_in == CEED_RESTRICTION_ORIENTED) { 1231 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_in, &orients_in)); 1232 } else if (rstr_type_in == CEED_RESTRICTION_CURL_ORIENTED) { 1233 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_in, &curl_orients_in)); 1234 } 1235 if (rstr_in != rstr_out) { 1236 if (rstr_type_out == CEED_RESTRICTION_ORIENTED) { 1237 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr_out, &orients_out)); 1238 } else if (rstr_type_out == CEED_RESTRICTION_CURL_ORIENTED) { 1239 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr_out, &curl_orients_out)); 1240 } 1241 } 1242 return CEED_ERROR_SUCCESS; 1243 } 1244 1245 //------------------------------------------------------------------------------ 1246 // Create operator 1247 //------------------------------------------------------------------------------ 1248 int CeedOperatorCreate_Cuda(CeedOperator op) { 1249 Ceed ceed; 1250 CeedOperator_Cuda *impl; 1251 1252 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 1253 CeedCallBackend(CeedCalloc(1, &impl)); 1254 CeedCallBackend(CeedOperatorSetData(op, impl)); 1255 1256 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Cuda)); 1257 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Cuda)); 1258 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddDiagonal", CeedOperatorLinearAssembleAddDiagonal_Cuda)); 1259 CeedCallBackend( 1260 CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleAddPointBlockDiagonal", CeedOperatorLinearAssembleAddPointBlockDiagonal_Cuda)); 1261 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleSingle", CeedSingleOperatorAssemble_Cuda)); 1262 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Cuda)); 1263 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Cuda)); 1264 return CEED_ERROR_SUCCESS; 1265 } 1266 1267 //------------------------------------------------------------------------------ 1268