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