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