1 // Copyright (c) 2017-2026, 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 <stdbool.h> 11 #include <stdint.h> 12 #include <string.h> 13 14 #include "ceed-opt.h" 15 16 //------------------------------------------------------------------------------ 17 // Setup Input/Output Fields 18 //------------------------------------------------------------------------------ 19 static int CeedOperatorSetupFields_Opt(CeedQFunction qf, CeedOperator op, bool is_input, bool *skip_rstr, bool *apply_add_basis, 20 const CeedInt block_size, CeedElemRestriction *block_rstr, CeedVector *e_vecs_full, CeedVector *e_vecs, 21 CeedVector *q_vecs, CeedInt start_e, CeedInt num_fields, CeedInt Q) { 22 Ceed ceed; 23 CeedSize e_size, q_size; 24 CeedInt num_comp, size, P; 25 CeedQFunctionField *qf_fields; 26 CeedOperatorField *op_fields; 27 28 { 29 Ceed ceed_parent; 30 31 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 32 CeedCallBackend(CeedGetParent(ceed, &ceed_parent)); 33 CeedCallBackend(CeedReferenceCopy(ceed_parent, &ceed)); 34 CeedCallBackend(CeedDestroy(&ceed_parent)); 35 } 36 if (is_input) { 37 CeedCallBackend(CeedOperatorGetFields(op, NULL, &op_fields, NULL, NULL)); 38 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_fields, NULL, NULL)); 39 } else { 40 CeedCallBackend(CeedOperatorGetFields(op, NULL, NULL, NULL, &op_fields)); 41 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, NULL, NULL, &qf_fields)); 42 } 43 44 // Loop over fields 45 for (CeedInt i = 0; i < num_fields; i++) { 46 CeedEvalMode eval_mode; 47 CeedBasis basis; 48 49 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_fields[i], &eval_mode)); 50 if (eval_mode != CEED_EVAL_WEIGHT) { 51 Ceed ceed_rstr; 52 CeedSize l_size; 53 CeedInt num_elem, elem_size, comp_stride; 54 CeedRestrictionType rstr_type; 55 CeedElemRestriction rstr; 56 57 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr)); 58 CeedCallBackend(CeedElemRestrictionGetCeed(rstr, &ceed_rstr)); 59 CeedCallBackend(CeedElemRestrictionGetNumElements(rstr, &num_elem)); 60 CeedCallBackend(CeedElemRestrictionGetElementSize(rstr, &elem_size)); 61 CeedCallBackend(CeedElemRestrictionGetLVectorSize(rstr, &l_size)); 62 CeedCallBackend(CeedElemRestrictionGetNumComponents(rstr, &num_comp)); 63 CeedCallBackend(CeedElemRestrictionGetCompStride(rstr, &comp_stride)); 64 65 CeedCallBackend(CeedElemRestrictionGetType(rstr, &rstr_type)); 66 switch (rstr_type) { 67 case CEED_RESTRICTION_STANDARD: { 68 const CeedInt *offsets = NULL; 69 70 CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets)); 71 CeedCallBackend(CeedElemRestrictionCreateBlocked(ceed_rstr, num_elem, elem_size, block_size, num_comp, comp_stride, l_size, CEED_MEM_HOST, 72 CEED_COPY_VALUES, offsets, &block_rstr[i + start_e])); 73 CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets)); 74 } break; 75 case CEED_RESTRICTION_ORIENTED: { 76 const bool *orients = NULL; 77 const CeedInt *offsets = NULL; 78 79 CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets)); 80 CeedCallBackend(CeedElemRestrictionGetOrientations(rstr, CEED_MEM_HOST, &orients)); 81 CeedCallBackend(CeedElemRestrictionCreateBlockedOriented(ceed_rstr, num_elem, elem_size, block_size, num_comp, comp_stride, l_size, 82 CEED_MEM_HOST, CEED_COPY_VALUES, offsets, orients, &block_rstr[i + start_e])); 83 CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets)); 84 CeedCallBackend(CeedElemRestrictionRestoreOrientations(rstr, &orients)); 85 } break; 86 case CEED_RESTRICTION_CURL_ORIENTED: { 87 const CeedInt8 *curl_orients = NULL; 88 const CeedInt *offsets = NULL; 89 90 CeedCallBackend(CeedElemRestrictionGetOffsets(rstr, CEED_MEM_HOST, &offsets)); 91 CeedCallBackend(CeedElemRestrictionGetCurlOrientations(rstr, CEED_MEM_HOST, &curl_orients)); 92 CeedCallBackend(CeedElemRestrictionCreateBlockedCurlOriented(ceed_rstr, num_elem, elem_size, block_size, num_comp, comp_stride, l_size, 93 CEED_MEM_HOST, CEED_COPY_VALUES, offsets, curl_orients, 94 &block_rstr[i + start_e])); 95 CeedCallBackend(CeedElemRestrictionRestoreOffsets(rstr, &offsets)); 96 CeedCallBackend(CeedElemRestrictionRestoreCurlOrientations(rstr, &curl_orients)); 97 } break; 98 case CEED_RESTRICTION_STRIDED: { 99 CeedInt strides[3]; 100 101 CeedCallBackend(CeedElemRestrictionGetStrides(rstr, strides)); 102 CeedCallBackend(CeedElemRestrictionCreateBlockedStrided(ceed_rstr, num_elem, elem_size, block_size, num_comp, l_size, strides, 103 &block_rstr[i + start_e])); 104 } break; 105 case CEED_RESTRICTION_POINTS: 106 // Empty case - won't occur 107 break; 108 } 109 CeedCallBackend(CeedDestroy(&ceed_rstr)); 110 CeedCallBackend(CeedElemRestrictionDestroy(&rstr)); 111 CeedCallBackend(CeedElemRestrictionCreateVector(block_rstr[i + start_e], NULL, &e_vecs_full[i + start_e])); 112 } 113 114 switch (eval_mode) { 115 case CEED_EVAL_NONE: 116 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 117 e_size = (CeedSize)Q * size * block_size; 118 CeedCallBackend(CeedVectorCreate(ceed, e_size, &e_vecs[i])); 119 q_size = (CeedSize)Q * size * block_size; 120 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 121 break; 122 case CEED_EVAL_INTERP: 123 case CEED_EVAL_GRAD: 124 case CEED_EVAL_DIV: 125 case CEED_EVAL_CURL: 126 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 127 CeedCallBackend(CeedQFunctionFieldGetSize(qf_fields[i], &size)); 128 CeedCallBackend(CeedBasisGetNumNodes(basis, &P)); 129 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 130 CeedCallBackend(CeedBasisDestroy(&basis)); 131 e_size = (CeedSize)P * num_comp * block_size; 132 CeedCallBackend(CeedVectorCreate(ceed, e_size, &e_vecs[i])); 133 q_size = (CeedSize)Q * size * block_size; 134 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 135 break; 136 case CEED_EVAL_WEIGHT: // Only on input fields 137 CeedCallBackend(CeedOperatorFieldGetBasis(op_fields[i], &basis)); 138 q_size = (CeedSize)Q * block_size; 139 CeedCallBackend(CeedVectorCreate(ceed, q_size, &q_vecs[i])); 140 CeedCallBackend(CeedBasisApply(basis, block_size, CEED_NOTRANSPOSE, CEED_EVAL_WEIGHT, CEED_VECTOR_NONE, q_vecs[i])); 141 CeedCallBackend(CeedBasisDestroy(&basis)); 142 break; 143 } 144 // Initialize E-vec arrays 145 if (e_vecs[i]) CeedCallBackend(CeedVectorSetValue(e_vecs[i], 0.0)); 146 } 147 // Drop duplicate restrictions 148 if (is_input) { 149 for (CeedInt i = 0; i < num_fields; i++) { 150 CeedVector vec_i; 151 CeedElemRestriction rstr_i; 152 153 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 154 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 155 for (CeedInt j = i + 1; j < num_fields; j++) { 156 CeedVector vec_j; 157 CeedElemRestriction rstr_j; 158 159 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 160 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 161 if (vec_i == vec_j && rstr_i == rstr_j) { 162 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j])); 163 CeedCallBackend(CeedVectorReferenceCopy(e_vecs_full[i + start_e], &e_vecs_full[j + start_e])); 164 skip_rstr[j] = true; 165 } 166 CeedCallBackend(CeedVectorDestroy(&vec_j)); 167 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j)); 168 } 169 CeedCallBackend(CeedVectorDestroy(&vec_i)); 170 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i)); 171 } 172 } else { 173 for (CeedInt i = num_fields - 1; i >= 0; i--) { 174 CeedVector vec_i; 175 CeedElemRestriction rstr_i; 176 177 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[i], &vec_i)); 178 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[i], &rstr_i)); 179 for (CeedInt j = i - 1; j >= 0; j--) { 180 CeedVector vec_j; 181 CeedElemRestriction rstr_j; 182 183 CeedCallBackend(CeedOperatorFieldGetVector(op_fields[j], &vec_j)); 184 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_fields[j], &rstr_j)); 185 if (vec_i == vec_j && rstr_i == rstr_j) { 186 CeedCallBackend(CeedVectorReferenceCopy(e_vecs[i], &e_vecs[j])); 187 CeedCallBackend(CeedVectorReferenceCopy(e_vecs_full[i + start_e], &e_vecs_full[j + start_e])); 188 skip_rstr[j] = true; 189 apply_add_basis[i] = true; 190 } 191 CeedCallBackend(CeedVectorDestroy(&vec_j)); 192 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_j)); 193 } 194 CeedCallBackend(CeedVectorDestroy(&vec_i)); 195 CeedCallBackend(CeedElemRestrictionDestroy(&rstr_i)); 196 } 197 } 198 CeedCallBackend(CeedDestroy(&ceed)); 199 return CEED_ERROR_SUCCESS; 200 } 201 202 //------------------------------------------------------------------------------ 203 // Setup Operator 204 //------------------------------------------------------------------------------ 205 static int CeedOperatorSetup_Opt(CeedOperator op) { 206 bool is_setup_done; 207 Ceed ceed; 208 Ceed_Opt *ceed_impl; 209 CeedInt Q, num_input_fields, num_output_fields; 210 CeedQFunctionField *qf_input_fields, *qf_output_fields; 211 CeedQFunction qf; 212 CeedOperatorField *op_input_fields, *op_output_fields; 213 CeedOperator_Opt *impl; 214 215 CeedCallBackend(CeedOperatorIsSetupDone(op, &is_setup_done)); 216 if (is_setup_done) return CEED_ERROR_SUCCESS; 217 218 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 219 CeedCallBackend(CeedGetData(ceed, &ceed_impl)); 220 CeedCallBackend(CeedDestroy(&ceed)); 221 CeedCallBackend(CeedOperatorGetData(op, &impl)); 222 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 223 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 224 CeedCallBackend(CeedQFunctionIsIdentity(qf, &impl->is_identity_qf)); 225 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 226 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 227 const CeedInt block_size = ceed_impl->block_size; 228 229 // Allocate 230 CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->block_rstr)); 231 CeedCallBackend(CeedCalloc(num_input_fields + num_output_fields, &impl->e_vecs_full)); 232 233 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_in)); 234 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->skip_rstr_out)); 235 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->apply_add_basis_out)); 236 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->input_states)); 237 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->e_vecs_in)); 238 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->e_vecs_out)); 239 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_in)); 240 CeedCallBackend(CeedCalloc(CEED_FIELD_MAX, &impl->q_vecs_out)); 241 242 impl->num_inputs = num_input_fields; 243 impl->num_outputs = num_output_fields; 244 245 // Set up infield and outfield pointer arrays 246 // Infields 247 CeedCallBackend(CeedOperatorSetupFields_Opt(qf, op, true, impl->skip_rstr_in, NULL, block_size, impl->block_rstr, impl->e_vecs_full, 248 impl->e_vecs_in, impl->q_vecs_in, 0, num_input_fields, Q)); 249 // Outfields 250 CeedCallBackend(CeedOperatorSetupFields_Opt(qf, op, false, impl->skip_rstr_out, impl->apply_add_basis_out, block_size, impl->block_rstr, 251 impl->e_vecs_full, impl->e_vecs_out, impl->q_vecs_out, num_input_fields, num_output_fields, Q)); 252 253 // Identity QFunctions 254 if (impl->is_identity_qf) { 255 CeedEvalMode in_mode, out_mode; 256 CeedQFunctionField *in_fields, *out_fields; 257 258 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &in_fields, NULL, &out_fields)); 259 CeedCallBackend(CeedQFunctionFieldGetEvalMode(in_fields[0], &in_mode)); 260 CeedCallBackend(CeedQFunctionFieldGetEvalMode(out_fields[0], &out_mode)); 261 262 if (in_mode == CEED_EVAL_NONE && out_mode == CEED_EVAL_NONE) { 263 impl->is_identity_rstr_op = true; 264 } else { 265 CeedCallBackend(CeedVectorReferenceCopy(impl->q_vecs_in[0], &impl->q_vecs_out[0])); 266 } 267 } 268 269 CeedCallBackend(CeedOperatorSetSetupDone(op)); 270 CeedCallBackend(CeedQFunctionDestroy(&qf)); 271 return CEED_ERROR_SUCCESS; 272 } 273 274 //------------------------------------------------------------------------------ 275 // Setup Input Fields 276 //------------------------------------------------------------------------------ 277 static inline int CeedOperatorSetupInputs_Opt(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 278 CeedVector in_vec, CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Opt *impl, 279 CeedRequest *request) { 280 for (CeedInt i = 0; i < num_input_fields; i++) { 281 CeedEvalMode eval_mode; 282 283 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 284 if (eval_mode == CEED_EVAL_WEIGHT) { // Skip 285 } else { 286 uint64_t state; 287 CeedVector vec; 288 289 // Get input vector 290 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 291 if (vec != CEED_VECTOR_ACTIVE) { 292 // Restrict 293 CeedCallBackend(CeedVectorGetState(vec, &state)); 294 if (state != impl->input_states[i] && impl->block_rstr[i] && !impl->skip_rstr_in[i]) { 295 CeedCallBackend(CeedElemRestrictionApply(impl->block_rstr[i], CEED_NOTRANSPOSE, vec, impl->e_vecs_full[i], request)); 296 } 297 impl->input_states[i] = state; 298 // Get evec 299 CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs_full[i], CEED_MEM_HOST, (const CeedScalar **)&e_data[i])); 300 } else { 301 // Set Qvec for CEED_EVAL_NONE 302 if (eval_mode == CEED_EVAL_NONE) { 303 CeedCallBackend(CeedVectorGetArrayRead(impl->e_vecs_in[i], CEED_MEM_HOST, (const CeedScalar **)&e_data[i])); 304 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_HOST, CEED_USE_POINTER, e_data[i])); 305 CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs_in[i], (const CeedScalar **)&e_data[i])); 306 } 307 } 308 CeedCallBackend(CeedVectorDestroy(&vec)); 309 } 310 } 311 return CEED_ERROR_SUCCESS; 312 } 313 314 //------------------------------------------------------------------------------ 315 // Input Basis Action 316 //------------------------------------------------------------------------------ 317 static inline int CeedOperatorInputBasis_Opt(CeedInt e, CeedInt Q, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 318 CeedInt num_input_fields, CeedInt block_size, CeedVector in_vec, bool skip_active, 319 CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Opt *impl, CeedRequest *request) { 320 for (CeedInt i = 0; i < num_input_fields; i++) { 321 bool is_active; 322 CeedInt elem_size, size, num_comp; 323 CeedEvalMode eval_mode; 324 CeedVector vec; 325 CeedElemRestriction elem_rstr; 326 CeedBasis basis; 327 328 // Skip active input 329 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 330 is_active = vec == CEED_VECTOR_ACTIVE; 331 CeedCallBackend(CeedVectorDestroy(&vec)); 332 if (skip_active && is_active) continue; 333 334 // Get elem_size, eval_mode, size 335 CeedCallBackend(CeedOperatorFieldGetElemRestriction(op_input_fields[i], &elem_rstr)); 336 CeedCallBackend(CeedElemRestrictionGetElementSize(elem_rstr, &elem_size)); 337 CeedCallBackend(CeedElemRestrictionDestroy(&elem_rstr)); 338 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 339 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &size)); 340 // Restrict block active input 341 if (is_active && impl->block_rstr[i]) { 342 CeedCallBackend(CeedElemRestrictionApplyBlock(impl->block_rstr[i], e / block_size, CEED_NOTRANSPOSE, in_vec, impl->e_vecs_in[i], request)); 343 } 344 // Basis action 345 switch (eval_mode) { 346 case CEED_EVAL_NONE: 347 if (!is_active) { 348 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_in[i], CEED_MEM_HOST, CEED_USE_POINTER, &e_data[i][(CeedSize)e * Q * size])); 349 } 350 break; 351 case CEED_EVAL_INTERP: 352 case CEED_EVAL_GRAD: 353 case CEED_EVAL_DIV: 354 case CEED_EVAL_CURL: 355 CeedCallBackend(CeedOperatorFieldGetBasis(op_input_fields[i], &basis)); 356 if (!is_active) { 357 CeedCallBackend(CeedBasisGetNumComponents(basis, &num_comp)); 358 CeedCallBackend(CeedVectorSetArray(impl->e_vecs_in[i], CEED_MEM_HOST, CEED_USE_POINTER, &e_data[i][(CeedSize)e * elem_size * num_comp])); 359 } 360 CeedCallBackend(CeedBasisApply(basis, block_size, CEED_NOTRANSPOSE, eval_mode, impl->e_vecs_in[i], impl->q_vecs_in[i])); 361 CeedCallBackend(CeedBasisDestroy(&basis)); 362 break; 363 case CEED_EVAL_WEIGHT: 364 break; // No action 365 } 366 } 367 return CEED_ERROR_SUCCESS; 368 } 369 370 //------------------------------------------------------------------------------ 371 // Output Basis Action 372 //------------------------------------------------------------------------------ 373 static inline int CeedOperatorOutputBasis_Opt(CeedInt e, CeedInt Q, CeedQFunctionField *qf_output_fields, CeedOperatorField *op_output_fields, 374 CeedInt block_size, CeedInt num_input_fields, CeedInt num_output_fields, bool *apply_add_basis, 375 bool *skip_rstr, CeedOperator op, CeedVector out_vec, CeedOperator_Opt *impl, CeedRequest *request) { 376 for (CeedInt i = 0; i < num_output_fields; i++) { 377 bool is_active; 378 CeedEvalMode eval_mode; 379 CeedVector vec; 380 CeedBasis basis; 381 382 // Get eval_mode 383 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 384 // Basis action 385 switch (eval_mode) { 386 case CEED_EVAL_NONE: 387 break; // No action 388 case CEED_EVAL_INTERP: 389 case CEED_EVAL_GRAD: 390 case CEED_EVAL_DIV: 391 case CEED_EVAL_CURL: 392 CeedCallBackend(CeedOperatorFieldGetBasis(op_output_fields[i], &basis)); 393 if (apply_add_basis[i]) { 394 CeedCallBackend(CeedBasisApplyAdd(basis, block_size, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs_out[i])); 395 } else { 396 CeedCallBackend(CeedBasisApply(basis, block_size, CEED_TRANSPOSE, eval_mode, impl->q_vecs_out[i], impl->e_vecs_out[i])); 397 } 398 CeedCallBackend(CeedBasisDestroy(&basis)); 399 break; 400 // LCOV_EXCL_START 401 case CEED_EVAL_WEIGHT: { 402 return CeedError(CeedOperatorReturnCeed(op), CEED_ERROR_BACKEND, "CEED_EVAL_WEIGHT cannot be an output evaluation mode"); 403 // LCOV_EXCL_STOP 404 } 405 } 406 // Restrict output block 407 if (skip_rstr[i]) continue; 408 // Get output vector 409 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 410 is_active = vec == CEED_VECTOR_ACTIVE; 411 if (is_active) vec = out_vec; 412 // Restrict 413 CeedCallBackend(CeedElemRestrictionApplyBlock(impl->block_rstr[i + impl->num_inputs], e / block_size, CEED_TRANSPOSE, impl->e_vecs_out[i], vec, 414 request)); 415 if (!is_active) CeedCallBackend(CeedVectorDestroy(&vec)); 416 } 417 return CEED_ERROR_SUCCESS; 418 } 419 420 //------------------------------------------------------------------------------ 421 // Restore Input Vectors 422 //------------------------------------------------------------------------------ 423 static inline int CeedOperatorRestoreInputs_Opt(CeedInt num_input_fields, CeedQFunctionField *qf_input_fields, CeedOperatorField *op_input_fields, 424 CeedScalar *e_data[2 * CEED_FIELD_MAX], CeedOperator_Opt *impl) { 425 for (CeedInt i = 0; i < num_input_fields; i++) { 426 CeedEvalMode eval_mode; 427 CeedVector vec; 428 429 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_input_fields[i], &eval_mode)); 430 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 431 if (eval_mode != CEED_EVAL_WEIGHT && vec != CEED_VECTOR_ACTIVE) { 432 CeedCallBackend(CeedVectorRestoreArrayRead(impl->e_vecs_full[i], (const CeedScalar **)&e_data[i])); 433 } 434 CeedCallBackend(CeedVectorDestroy(&vec)); 435 } 436 return CEED_ERROR_SUCCESS; 437 } 438 439 //------------------------------------------------------------------------------ 440 // Operator Apply 441 //------------------------------------------------------------------------------ 442 static int CeedOperatorApplyAdd_Opt(CeedOperator op, CeedVector in_vec, CeedVector out_vec, CeedRequest *request) { 443 Ceed ceed; 444 Ceed_Opt *ceed_impl; 445 CeedInt Q, num_input_fields, num_output_fields, num_elem; 446 CeedEvalMode eval_mode; 447 CeedScalar *e_data[2 * CEED_FIELD_MAX] = {0}; 448 CeedQFunctionField *qf_input_fields, *qf_output_fields; 449 CeedQFunction qf; 450 CeedOperatorField *op_input_fields, *op_output_fields; 451 CeedOperator_Opt *impl; 452 453 // Setup 454 CeedCallBackend(CeedOperatorSetup_Opt(op)); 455 456 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 457 CeedCallBackend(CeedGetData(ceed, &ceed_impl)); 458 CeedCallBackend(CeedDestroy(&ceed)); 459 CeedCallBackend(CeedOperatorGetData(op, &impl)); 460 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 461 const CeedInt block_size = ceed_impl->block_size; 462 const CeedInt num_blocks = (num_elem / block_size) + !!(num_elem % block_size); 463 464 // Restriction only operator 465 if (impl->is_identity_rstr_op) { 466 for (CeedInt b = 0; b < num_blocks; b++) { 467 CeedCallBackend(CeedElemRestrictionApplyBlock(impl->block_rstr[0], b, CEED_NOTRANSPOSE, in_vec, impl->e_vecs_in[0], request)); 468 CeedCallBackend(CeedElemRestrictionApplyBlock(impl->block_rstr[1], b, CEED_TRANSPOSE, impl->e_vecs_in[0], out_vec, request)); 469 } 470 return CEED_ERROR_SUCCESS; 471 } 472 473 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 474 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 475 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 476 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 477 478 // Input Evecs and Restriction 479 CeedCallBackend(CeedOperatorSetupInputs_Opt(num_input_fields, qf_input_fields, op_input_fields, in_vec, e_data, impl, request)); 480 481 // Output Lvecs, Evecs, and Qvecs 482 for (CeedInt i = 0; i < num_output_fields; i++) { 483 // Set Qvec if needed 484 CeedCallBackend(CeedQFunctionFieldGetEvalMode(qf_output_fields[i], &eval_mode)); 485 if (eval_mode == CEED_EVAL_NONE) { 486 // Set qvec to single block evec 487 CeedCallBackend(CeedVectorGetArrayWrite(impl->e_vecs_out[i], CEED_MEM_HOST, &e_data[i + num_input_fields])); 488 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[i], CEED_MEM_HOST, CEED_USE_POINTER, e_data[i + num_input_fields])); 489 CeedCallBackend(CeedVectorRestoreArray(impl->e_vecs_out[i], &e_data[i + num_input_fields])); 490 } 491 } 492 493 // Loop through elements 494 for (CeedInt e = 0; e < num_blocks * block_size; e += block_size) { 495 // Input basis apply 496 CeedCallBackend(CeedOperatorInputBasis_Opt(e, Q, qf_input_fields, op_input_fields, num_input_fields, block_size, in_vec, false, e_data, impl, 497 request)); 498 499 // Q function 500 if (!impl->is_identity_qf) { 501 CeedCallBackend(CeedQFunctionApply(qf, Q * block_size, impl->q_vecs_in, impl->q_vecs_out)); 502 } 503 504 // Output basis apply and restriction 505 CeedCallBackend(CeedOperatorOutputBasis_Opt(e, Q, qf_output_fields, op_output_fields, block_size, num_input_fields, num_output_fields, 506 impl->apply_add_basis_out, impl->skip_rstr_out, op, out_vec, impl, request)); 507 } 508 509 // Restore input arrays 510 CeedCallBackend(CeedOperatorRestoreInputs_Opt(num_input_fields, qf_input_fields, op_input_fields, e_data, impl)); 511 CeedCallBackend(CeedQFunctionDestroy(&qf)); 512 return CEED_ERROR_SUCCESS; 513 } 514 515 //------------------------------------------------------------------------------ 516 // Core code for linear QFunction assembly 517 //------------------------------------------------------------------------------ 518 static inline int CeedOperatorLinearAssembleQFunctionCore_Opt(CeedOperator op, bool build_objects, CeedVector *assembled, CeedElemRestriction *rstr, 519 CeedRequest *request) { 520 Ceed ceed; 521 Ceed_Opt *ceed_impl; 522 CeedInt qf_size_in, qf_size_out, Q, num_input_fields, num_output_fields, num_elem; 523 CeedScalar *l_vec_array, *e_data[2 * CEED_FIELD_MAX] = {0}; 524 CeedQFunctionField *qf_input_fields, *qf_output_fields; 525 CeedQFunction qf; 526 CeedOperatorField *op_input_fields, *op_output_fields; 527 CeedOperator_Opt *impl; 528 529 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 530 CeedCallBackend(CeedGetData(ceed, &ceed_impl)); 531 CeedCallBackend(CeedOperatorGetData(op, &impl)); 532 qf_size_in = impl->qf_size_in; 533 qf_size_out = impl->qf_size_out; 534 535 CeedCallBackend(CeedOperatorGetNumElements(op, &num_elem)); 536 CeedCallBackend(CeedOperatorGetNumQuadraturePoints(op, &Q)); 537 CeedCallBackend(CeedOperatorGetQFunction(op, &qf)); 538 CeedCallBackend(CeedOperatorGetFields(op, &num_input_fields, &op_input_fields, &num_output_fields, &op_output_fields)); 539 CeedCallBackend(CeedQFunctionGetFields(qf, NULL, &qf_input_fields, NULL, &qf_output_fields)); 540 const CeedInt block_size = ceed_impl->block_size; 541 const CeedInt num_blocks = (num_elem / block_size) + !!(num_elem % block_size); 542 CeedVector l_vec = impl->qf_l_vec; 543 CeedElemRestriction block_rstr = impl->qf_block_rstr; 544 545 // Setup 546 CeedCallBackend(CeedOperatorSetup_Opt(op)); 547 548 // Check for restriction only operator 549 CeedCheck(!impl->is_identity_rstr_op, ceed, CEED_ERROR_BACKEND, "Assembling restriction only operators is not supported"); 550 551 // Input Evecs and Restriction 552 CeedCallBackend(CeedOperatorSetupInputs_Opt(num_input_fields, qf_input_fields, op_input_fields, NULL, e_data, impl, request)); 553 554 // Count number of active input fields 555 if (qf_size_in == 0) { 556 for (CeedInt i = 0; i < num_input_fields; i++) { 557 CeedInt field_size; 558 CeedVector vec; 559 560 // Check if active input 561 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 562 if (vec == CEED_VECTOR_ACTIVE) { 563 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &field_size)); 564 CeedCallBackend(CeedVectorSetValue(impl->q_vecs_in[i], 0.0)); 565 qf_size_in += field_size; 566 } 567 CeedCallBackend(CeedVectorDestroy(&vec)); 568 } 569 CeedCheck(qf_size_in > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 570 impl->qf_size_in = qf_size_in; 571 } 572 573 // Count number of active output fields 574 if (qf_size_out == 0) { 575 for (CeedInt i = 0; i < num_output_fields; i++) { 576 CeedInt field_size; 577 CeedVector vec; 578 579 // Check if active output 580 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[i], &vec)); 581 if (vec == CEED_VECTOR_ACTIVE) { 582 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[i], &field_size)); 583 qf_size_out += field_size; 584 } 585 CeedCallBackend(CeedVectorDestroy(&vec)); 586 } 587 CeedCheck(qf_size_out > 0, ceed, CEED_ERROR_BACKEND, "Cannot assemble QFunction without active inputs and outputs"); 588 impl->qf_size_out = qf_size_out; 589 } 590 591 // Setup l_vec 592 if (!l_vec) { 593 const CeedSize l_size = (CeedSize)block_size * Q * qf_size_in * qf_size_out; 594 595 CeedCallBackend(CeedVectorCreate(ceed, l_size, &l_vec)); 596 CeedCallBackend(CeedVectorSetValue(l_vec, 0.0)); 597 impl->qf_l_vec = l_vec; 598 } 599 600 // Output blocked restriction 601 if (!block_rstr) { 602 CeedInt strides[3] = {1, Q, qf_size_in * qf_size_out * Q}; 603 604 CeedCallBackend(CeedElemRestrictionCreateBlockedStrided(ceed, num_elem, Q, block_size, qf_size_in * qf_size_out, 605 qf_size_in * qf_size_out * num_elem * Q, strides, &block_rstr)); 606 impl->qf_block_rstr = block_rstr; 607 } 608 609 // Build objects if needed 610 if (build_objects) { 611 const CeedSize l_size = (CeedSize)num_elem * Q * qf_size_in * qf_size_out; 612 CeedInt strides[3] = {1, Q, qf_size_in * qf_size_out * Q}; 613 614 // Create output restriction 615 CeedCallBackend(CeedElemRestrictionCreateStrided(ceed, num_elem, Q, qf_size_in * qf_size_out, 616 (CeedSize)qf_size_in * (CeedSize)qf_size_out * (CeedSize)num_elem * (CeedSize)Q, strides, rstr)); 617 // Create assembled vector 618 CeedCallBackend(CeedVectorCreate(ceed, l_size, assembled)); 619 } 620 621 // Loop through elements 622 CeedCallBackend(CeedVectorSetValue(*assembled, 0.0)); 623 for (CeedInt e = 0; e < num_blocks * block_size; e += block_size) { 624 CeedCallBackend(CeedVectorGetArray(l_vec, CEED_MEM_HOST, &l_vec_array)); 625 626 // Input basis apply 627 CeedCallBackend(CeedOperatorInputBasis_Opt(e, Q, qf_input_fields, op_input_fields, num_input_fields, block_size, NULL, true, e_data, impl, 628 request)); 629 630 // Assemble QFunction 631 for (CeedInt i = 0; i < num_input_fields; i++) { 632 bool is_active; 633 CeedInt field_size; 634 CeedVector vec; 635 636 // Check if active input 637 CeedCallBackend(CeedOperatorFieldGetVector(op_input_fields[i], &vec)); 638 is_active = vec == CEED_VECTOR_ACTIVE; 639 CeedCallBackend(CeedVectorDestroy(&vec)); 640 if (!is_active) continue; 641 CeedCallBackend(CeedQFunctionFieldGetSize(qf_input_fields[i], &field_size)); 642 for (CeedInt field = 0; field < field_size; field++) { 643 // Set current portion of input to 1.0 644 { 645 CeedScalar *array; 646 647 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_HOST, &array)); 648 for (CeedInt j = 0; j < Q * block_size; j++) array[field * Q * block_size + j] = 1.0; 649 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &array)); 650 } 651 652 if (!impl->is_identity_qf) { 653 // Set Outputs 654 for (CeedInt out = 0; out < num_output_fields; out++) { 655 CeedVector vec; 656 657 // Check if active output 658 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 659 if (vec == CEED_VECTOR_ACTIVE) { 660 CeedInt field_size; 661 662 CeedCallBackend(CeedVectorSetArray(impl->q_vecs_out[out], CEED_MEM_HOST, CEED_USE_POINTER, l_vec_array)); 663 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[out], &field_size)); 664 l_vec_array += field_size * Q * block_size; // Advance the pointer by the size of the output 665 } 666 CeedCallBackend(CeedVectorDestroy(&vec)); 667 } 668 // Apply QFunction 669 CeedCallBackend(CeedQFunctionApply(qf, Q * block_size, impl->q_vecs_in, impl->q_vecs_out)); 670 } else { 671 CeedInt field_size; 672 const CeedScalar *array; 673 674 // Copy Identity Outputs 675 CeedCallBackend(CeedQFunctionFieldGetSize(qf_output_fields[0], &field_size)); 676 CeedCallBackend(CeedVectorGetArrayRead(impl->q_vecs_out[0], CEED_MEM_HOST, &array)); 677 for (CeedInt j = 0; j < field_size * Q * block_size; j++) l_vec_array[j] = array[j]; 678 CeedCallBackend(CeedVectorRestoreArrayRead(impl->q_vecs_out[0], &array)); 679 l_vec_array += field_size * Q * block_size; 680 } 681 // Reset input to 0.0 682 { 683 CeedScalar *array; 684 685 CeedCallBackend(CeedVectorGetArray(impl->q_vecs_in[i], CEED_MEM_HOST, &array)); 686 for (CeedInt j = 0; j < Q * block_size; j++) array[field * Q * block_size + j] = 0.0; 687 CeedCallBackend(CeedVectorRestoreArray(impl->q_vecs_in[i], &array)); 688 } 689 } 690 } 691 692 // Un-set output Qvecs to prevent accidental overwrite of Assembled 693 if (!impl->is_identity_qf) { 694 for (CeedInt out = 0; out < num_output_fields; out++) { 695 CeedVector vec; 696 697 // Check if active output 698 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 699 if (vec == CEED_VECTOR_ACTIVE && num_elem > 0) { 700 CeedCallBackend(CeedVectorTakeArray(impl->q_vecs_out[out], CEED_MEM_HOST, NULL)); 701 } 702 CeedCallBackend(CeedVectorDestroy(&vec)); 703 } 704 } 705 706 // Assemble into assembled vector 707 CeedCallBackend(CeedVectorRestoreArray(l_vec, &l_vec_array)); 708 CeedCallBackend(CeedElemRestrictionApplyBlock(block_rstr, e / block_size, CEED_TRANSPOSE, l_vec, *assembled, request)); 709 } 710 711 // Reset output Qvecs 712 for (CeedInt out = 0; out < num_output_fields; out++) { 713 CeedVector vec; 714 715 // Initialize array if active output 716 CeedCallBackend(CeedOperatorFieldGetVector(op_output_fields[out], &vec)); 717 if (vec == CEED_VECTOR_ACTIVE) CeedCallBackend(CeedVectorSetValue(impl->q_vecs_out[out], 0.0)); 718 CeedCallBackend(CeedVectorDestroy(&vec)); 719 } 720 721 // Restore input arrays 722 CeedCallBackend(CeedOperatorRestoreInputs_Opt(num_input_fields, qf_input_fields, op_input_fields, e_data, impl)); 723 CeedCallBackend(CeedDestroy(&ceed)); 724 CeedCallBackend(CeedQFunctionDestroy(&qf)); 725 return CEED_ERROR_SUCCESS; 726 } 727 728 //------------------------------------------------------------------------------ 729 // Assemble Linear QFunction 730 //------------------------------------------------------------------------------ 731 static int CeedOperatorLinearAssembleQFunction_Opt(CeedOperator op, CeedVector *assembled, CeedElemRestriction *rstr, CeedRequest *request) { 732 return CeedOperatorLinearAssembleQFunctionCore_Opt(op, true, assembled, rstr, request); 733 } 734 735 //------------------------------------------------------------------------------ 736 // Update Assembled Linear QFunction 737 //------------------------------------------------------------------------------ 738 static int CeedOperatorLinearAssembleQFunctionUpdate_Opt(CeedOperator op, CeedVector assembled, CeedElemRestriction rstr, CeedRequest *request) { 739 return CeedOperatorLinearAssembleQFunctionCore_Opt(op, false, &assembled, &rstr, request); 740 } 741 742 //------------------------------------------------------------------------------ 743 // Operator Destroy 744 //------------------------------------------------------------------------------ 745 static int CeedOperatorDestroy_Opt(CeedOperator op) { 746 CeedOperator_Opt *impl; 747 748 CeedCallBackend(CeedOperatorGetData(op, &impl)); 749 for (CeedInt i = 0; i < impl->num_inputs + impl->num_outputs; i++) { 750 CeedCallBackend(CeedElemRestrictionDestroy(&impl->block_rstr[i])); 751 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_full[i])); 752 } 753 CeedCallBackend(CeedFree(&impl->block_rstr)); 754 CeedCallBackend(CeedFree(&impl->e_vecs_full)); 755 CeedCallBackend(CeedFree(&impl->input_states)); 756 CeedCallBackend(CeedFree(&impl->skip_rstr_in)); 757 CeedCallBackend(CeedFree(&impl->skip_rstr_out)); 758 CeedCallBackend(CeedFree(&impl->apply_add_basis_out)); 759 760 for (CeedInt i = 0; i < impl->num_inputs; i++) { 761 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_in[i])); 762 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_in[i])); 763 } 764 CeedCallBackend(CeedFree(&impl->e_vecs_in)); 765 CeedCallBackend(CeedFree(&impl->q_vecs_in)); 766 767 for (CeedInt i = 0; i < impl->num_outputs; i++) { 768 CeedCallBackend(CeedVectorDestroy(&impl->e_vecs_out[i])); 769 CeedCallBackend(CeedVectorDestroy(&impl->q_vecs_out[i])); 770 } 771 CeedCallBackend(CeedFree(&impl->e_vecs_out)); 772 CeedCallBackend(CeedFree(&impl->q_vecs_out)); 773 774 // QFunction assembly data 775 CeedCallBackend(CeedVectorDestroy(&impl->qf_l_vec)); 776 CeedCallBackend(CeedElemRestrictionDestroy(&impl->qf_block_rstr)); 777 778 CeedCallBackend(CeedFree(&impl)); 779 return CEED_ERROR_SUCCESS; 780 } 781 782 //------------------------------------------------------------------------------ 783 // Operator Create 784 //------------------------------------------------------------------------------ 785 int CeedOperatorCreate_Opt(CeedOperator op) { 786 Ceed ceed; 787 Ceed_Opt *ceed_impl; 788 CeedOperator_Opt *impl; 789 790 CeedCallBackend(CeedOperatorGetCeed(op, &ceed)); 791 CeedCallBackend(CeedGetData(ceed, &ceed_impl)); 792 const CeedInt block_size = ceed_impl->block_size; 793 794 CeedCallBackend(CeedCalloc(1, &impl)); 795 CeedCallBackend(CeedOperatorSetData(op, impl)); 796 797 CeedCheck(block_size == 1 || block_size == 8, ceed, CEED_ERROR_BACKEND, "Opt backend cannot use blocksize: %" CeedInt_FMT, block_size); 798 799 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunction", CeedOperatorLinearAssembleQFunction_Opt)); 800 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "LinearAssembleQFunctionUpdate", CeedOperatorLinearAssembleQFunctionUpdate_Opt)); 801 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "ApplyAdd", CeedOperatorApplyAdd_Opt)); 802 CeedCallBackend(CeedSetBackendFunction(ceed, "Operator", op, "Destroy", CeedOperatorDestroy_Opt)); 803 CeedCallBackend(CeedDestroy(&ceed)); 804 return CEED_ERROR_SUCCESS; 805 } 806 807 //------------------------------------------------------------------------------ 808