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