1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 362b7942eSJames Wright 462b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h" 562b7942eSJames Wright 662b7942eSJames Wright #include <petscdmplex.h> 762b7942eSJames Wright 8149fb536SJames Wright #include <navierstokes.h> 94c07ec22SJames Wright #include <sgs_model_torch.h> 1062b7942eSJames Wright 1182baf964SJames Wright typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 1282baf964SJames Wright typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 1382baf964SJames Wright typedef struct { 1482baf964SJames Wright DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 1582baf964SJames Wright PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 1682baf964SJames Wright OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 1782baf964SJames Wright CeedVector sgs_nodal_ceed, grad_velo_ceed; 1882baf964SJames Wright SgsDDNodalStressEval sgs_nodal_eval; 1982baf964SJames Wright SgsDDNodalStressInference sgs_nodal_inference; 2082baf964SJames Wright void *sgs_nodal_inference_ctx; 2182baf964SJames Wright PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 2282baf964SJames Wright } *SgsDDData; 2382baf964SJames Wright 2482baf964SJames Wright // @brief Destroy `SgsDDData` object 2582baf964SJames Wright static PetscErrorCode SgsDDDataDestroy(SgsDDData *sgs_dd_data) { 2682baf964SJames Wright SgsDDData sgs_dd_data_ = *sgs_dd_data; 27*14bd2a07SJames Wright 2882baf964SJames Wright PetscFunctionBeginUser; 2982baf964SJames Wright if (!sgs_dd_data_) PetscFunctionReturn(PETSC_SUCCESS); 3082baf964SJames Wright Ceed ceed = sgs_dd_data_->op_sgs_apply_ctx->ceed; 3182baf964SJames Wright 3282baf964SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->sgs_nodal_ceed)); 3382baf964SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->grad_velo_ceed)); 3482baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_evaluation_ctx)); 3582baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_sgs_apply_ctx)); 3682baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_inputs_ctx)); 3782baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_outputs_ctx)); 3882baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_sgs)); 3982baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_inputs)); 4082baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_outputs)); 4182baf964SJames Wright if (sgs_dd_data_->sgs_nodal_inference_ctx) PetscCall(sgs_dd_data_->sgs_nodal_inference_ctx_destroy(sgs_dd_data_->sgs_nodal_inference_ctx)); 4282baf964SJames Wright PetscCall(PetscFree(sgs_dd_data_)); 4382baf964SJames Wright *sgs_dd_data = NULL; 4482baf964SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4582baf964SJames Wright } 4682baf964SJames Wright 47c38c977aSJames Wright typedef struct { 48c38c977aSJames Wright CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 49c38c977aSJames Wright CeedVector grid_aniso_ceed; 5040816385SJames Wright CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx; 514c07ec22SJames Wright SGSModelDDImplementation sgs_dd_model_implementation; 52ad494f68SJames Wright } *SgsDDSetupData; 53c38c977aSJames Wright 548340219bSJames Wright #define GRAD_VELO_PROJ_KEY "Gradient of Velocity Projection" 5582baf964SJames Wright #define SGS_DD_DATA_KEY "SGS Data Driven Data" 568340219bSJames Wright 57ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) { 58b4c37c5cSJames Wright Ceed ceed; 59ad494f68SJames Wright 60c38c977aSJames Wright PetscFunctionBeginUser; 61b4c37c5cSJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 62ad494f68SJames Wright 63b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 64b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 65b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 66b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 6740816385SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx)); 68c38c977aSJames Wright PetscCall(PetscFree(sgs_dd_setup_data)); 69519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 70d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 71c38c977aSJames Wright } 72c38c977aSJames Wright 73ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes 74ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 75ee1455b7SJames Wright PetscSection section; 76ee1455b7SJames Wright 77ee1455b7SJames Wright PetscFunctionBeginUser; 78ee1455b7SJames Wright *num_components = 6; 79ee1455b7SJames Wright 80ee1455b7SJames Wright PetscCall(DMClone(dm_source, dm_sgs)); 810dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE)); 82ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 83ee1455b7SJames Wright 84da4ca0cfSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 85ee1455b7SJames Wright 86ee1455b7SJames Wright PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 87ee1455b7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 88ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 89ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 90ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 91ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 92ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 93ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 94d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 95ee1455b7SJames Wright }; 96ee1455b7SJames Wright 97ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method 980c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 9982baf964SJames Wright SgsDDData sgs_dd_data; 100cceb3143SJames Wright PetscMemType q_mem_type; 101cceb3143SJames Wright 102cceb3143SJames Wright PetscFunctionBeginUser; 1030c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 1040c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 105cceb3143SJames Wright 106cceb3143SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 107cceb3143SJames Wright 1080c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 109cceb3143SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 110cceb3143SJames Wright } 111cceb3143SJames Wright 112b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator 113e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 11482baf964SJames Wright SgsDDData sgs_dd_data; 1155930f037SJames Wright CeedQFunction qf_sgs_dd_nodal; 1165930f037SJames Wright CeedOperator op_sgs_dd_nodal; 11715c18037SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 118defe8520SJames Wright PetscInt dim; 1195930f037SJames Wright CeedVector inv_multiplicity; 120ee1455b7SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 12115c18037SJames Wright DMLabel domain_label = NULL; 12215c18037SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 1238340219bSJames Wright NodalProjectionData grad_velo_proj; 124ee1455b7SJames Wright 125ee1455b7SJames Wright PetscFunctionBeginUser; 1260c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 127e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 128e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 129b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 1300c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 1310c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 132ee1455b7SJames Wright 133ee1455b7SJames Wright { // Get velocity gradient information 134ee1455b7SJames Wright CeedOperatorField op_field; 1358340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 137b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 138ee1455b7SJames Wright } 13915c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 140b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 141ee1455b7SJames Wright 1425930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 1435930f037SJames Wright &inv_multiplicity)); 144ee1455b7SJames Wright 145ee1455b7SJames Wright // -- Create operator for SGS DD model nodal evaluation 1460c373b74SJames Wright switch (honee->phys->state_var) { 147ee1455b7SJames Wright case STATEVAR_PRIMITIVE: 148ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal)); 149ee1455b7SJames Wright break; 150ee1455b7SJames Wright case STATEVAR_CONSERVATIVE: 151ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal)); 152ee1455b7SJames Wright break; 1539b103f75SJames Wright case STATEVAR_ENTROPY: 1549b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal)); 1559b103f75SJames Wright break; 156ee1455b7SJames Wright } 157ee1455b7SJames Wright 158ee1455b7SJames Wright // Mesh/geometry order and solution basis order may differ, therefore must interpolate 159ee1455b7SJames Wright CeedBasis basis_x_to_q; 160e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q)); 161ee1455b7SJames Wright 162b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 163b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 164b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 165b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 166b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 167b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 168b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 169ee1455b7SJames Wright 170b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 171e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 172e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord)); 17358e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 17458e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 175b4c37c5cSJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 17658e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 17758e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 178ee1455b7SJames Wright 1798340219bSJames Wright PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, NULL, sgs_dd_data->sgs_nodal_ceed, NULL, NULL, 1808340219bSJames Wright &sgs_dd_data->op_nodal_evaluation_ctx)); 181ee1455b7SJames Wright 182ee1455b7SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 183ad494f68SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Fused; 184ee1455b7SJames Wright 185b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 186b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 187b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 188fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 189b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 190b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 191d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 192ee1455b7SJames Wright } 193ee1455b7SJames Wright 1944c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation 1954c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 1964c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 197b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 198b87d60b3SJames Wright void **ctx) { 199b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inference; 200b87d60b3SJames Wright CeedOperator op_sgs_dd_inference; 201b87d60b3SJames Wright OperatorApplyContext *op_context = (OperatorApplyContext *)ctx; 202b87d60b3SJames Wright 203b87d60b3SJames Wright PetscFunctionBeginUser; 204b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc, 205b87d60b3SJames Wright &qf_sgs_dd_inference)); 206b87d60b3SJames Wright 207b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx)); 208b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 209b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE)); 210b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 211b87d60b3SJames Wright 212b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference)); 213b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 214b87d60b3SJames Wright PetscCallCeed(ceed, 215b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 216b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 217b87d60b3SJames Wright 218b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL, 219b87d60b3SJames Wright op_context)); 220b87d60b3SJames Wright sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy; 221b87d60b3SJames Wright 222b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference)); 223b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference)); 224b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 225b87d60b3SJames Wright } 226b87d60b3SJames Wright 2274c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation 2284c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 229b87d60b3SJames Wright OperatorApplyContext op_context = *(OperatorApplyContext *)ctx; 230b87d60b3SJames Wright 231b87d60b3SJames Wright PetscFunctionBeginUser; 232ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 233ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 234b40a7e63SJames Wright PetscCall(PetscLogGpuTimeBegin()); 235b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context)); 236b40a7e63SJames Wright PetscCall(PetscLogGpuTimeEnd()); 237ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 238ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 239b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 240b87d60b3SJames Wright } 241b87d60b3SJames Wright 2424c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch 2434c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 2444c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 2454c07ec22SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 2464c07ec22SJames Wright void **ctx) { 2474c07ec22SJames Wright const char *ceed_resource; 2484c07ec22SJames Wright char model_path[PETSC_MAX_PATH_LEN] = ""; 2494c07ec22SJames Wright TorchDeviceType model_device_type; 2504c07ec22SJames Wright 2514c07ec22SJames Wright PetscFunctionBeginUser; 2524c07ec22SJames Wright PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource)); 2534c07ec22SJames Wright if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA; 2544c07ec22SJames Wright else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP; 2557ffa0ff8SJames Wright // On-device XPU is not working reliably currently, default to CPU inference evaluation 2567ffa0ff8SJames Wright // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU; 2574c07ec22SJames Wright else model_device_type = TORCH_DEVICE_CPU; 2584c07ec22SJames Wright PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL)); 2594c07ec22SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL)); 2604c07ec22SJames Wright 2614c07ec22SJames Wright PetscCall(LoadModel_Torch(model_path, model_device_type)); 2624c07ec22SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2634c07ec22SJames Wright } 2644c07ec22SJames Wright 2654c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch 2664c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 2674c07ec22SJames Wright static PetscBool run_through = PETSC_FALSE; 268*14bd2a07SJames Wright 2694c07ec22SJames Wright PetscFunctionBeginUser; 2704c07ec22SJames Wright if (!run_through) { 2714c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view")); 2724c07ec22SJames Wright } 2734c07ec22SJames Wright PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc)); 2744c07ec22SJames Wright if (!run_through) { 2754c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view")); 2764c07ec22SJames Wright run_through = PETSC_TRUE; 2774c07ec22SJames Wright } 2784c07ec22SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2794c07ec22SJames Wright } 2804c07ec22SJames Wright 281b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method 2820c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 28382baf964SJames Wright SgsDDData sgs_dd_data; 284b87d60b3SJames Wright PetscMemType q_mem_type; 285b87d60b3SJames Wright Vec DD_Inputs_loc, DD_Outputs_loc; 286b87d60b3SJames Wright 287b87d60b3SJames Wright PetscFunctionBeginUser; 2880c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 289b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 290b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 2910c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 292b87d60b3SJames Wright 293b87d60b3SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx)); 294b87d60b3SJames Wright PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx)); 295b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx)); 296b87d60b3SJames Wright 2970c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 298b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 299b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 300b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 301b87d60b3SJames Wright } 302b87d60b3SJames Wright 303b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators 304e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 30582baf964SJames Wright SgsDDData sgs_dd_data; 306b87d60b3SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1; 307b87d60b3SJames Wright PetscInt dim; 308b87d60b3SJames Wright CeedVector inv_multiplicity, eigvec; 3098340219bSJames Wright NodalProjectionData grad_velo_proj; 310b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs, 311b87d60b3SJames Wright elem_restr_dd_outputs; 312b87d60b3SJames Wright DMLabel domain_label = NULL; 313b87d60b3SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 314b87d60b3SJames Wright 315b87d60b3SJames Wright PetscFunctionBeginUser; 3160c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 317b87d60b3SJames Wright { // Create DMs for data-driven input and output values 318b87d60b3SJames Wright PetscSection section; 319b87d60b3SJames Wright PetscInt degree, q_extra; 320b87d60b3SJames Wright { // Get degree and number of quadrature points from dm_sgs 321b87d60b3SJames Wright PetscFE fe; 322b87d60b3SJames Wright PetscSpace basis; 323b87d60b3SJames Wright PetscQuadrature quadrature; 324b87d60b3SJames Wright PetscInt num_qpnts; 325b87d60b3SJames Wright PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe)); 326b87d60b3SJames Wright PetscCall(PetscFEGetBasisSpace(fe, &basis)); 327b87d60b3SJames Wright PetscCall(PetscSpaceGetDegree(basis, °ree, NULL)); 328b87d60b3SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 329b87d60b3SJames Wright PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts)); 330b87d60b3SJames Wright q_extra = degree - num_qpnts; 331b87d60b3SJames Wright } 332b87d60b3SJames Wright 333b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs)); 3340dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE)); 335b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs")); 336b87d60b3SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_inputs, sgs_dd_data->dm_dd_inputs)); 337b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, §ion)); 338b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 339b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) { 340b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 341b87d60b3SJames Wright 342b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1)); 343b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 344b87d60b3SJames Wright } 345b87d60b3SJames Wright 346b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs)); 3470dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE)); 348b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs")); 349b87d60b3SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_outputs, sgs_dd_data->dm_dd_outputs)); 350b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, §ion)); 351b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 352b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) { 353b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 354b87d60b3SJames Wright 355b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1)); 356b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 357b87d60b3SJames Wright } 358b87d60b3SJames Wright } 359b87d60b3SJames Wright 3600c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 361e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 362e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 363b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 3640c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 365b87d60b3SJames Wright 366b87d60b3SJames Wright { // Get velocity gradient information 367b87d60b3SJames Wright CeedOperatorField op_field; 3688340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 369b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 370b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 371b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL)); 372b87d60b3SJames Wright } 373b87d60b3SJames Wright 374b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 375b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 376b87d60b3SJames Wright PetscCall( 377b87d60b3SJames Wright DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec)); 378b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL)); 379b87d60b3SJames Wright 380b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs)); 381b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs)); 382b87d60b3SJames Wright 3835930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 3845930f037SJames Wright &inv_multiplicity)); 385b87d60b3SJames Wright 386b87d60b3SJames Wright { // Create operator for data-driven input evaluation 387b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inputs; 388b87d60b3SJames Wright CeedOperator op_sgs_dd_inputs; 389b87d60b3SJames Wright 3900c373b74SJames Wright switch (honee->phys->state_var) { 391b87d60b3SJames Wright case STATEVAR_PRIMITIVE: 392b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim, 393b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs)); 394b87d60b3SJames Wright break; 395b87d60b3SJames Wright case STATEVAR_CONSERVATIVE: 396b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv, 397b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs)); 398b87d60b3SJames Wright break; 3999b103f75SJames Wright case STATEVAR_ENTROPY: 4009b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy, 4019b103f75SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs)); 4029b103f75SJames Wright break; 403b87d60b3SJames Wright } 404b87d60b3SJames Wright 405b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx)); 406b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE)); 407b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 408b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 409b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 410b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 411b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 412b87d60b3SJames Wright 413b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs)); 414e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 415b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 416b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 417b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 418b87d60b3SJames Wright PetscCallCeed(ceed, 419b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 420b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 421b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 422b87d60b3SJames Wright 4238340219bSJames Wright PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL, 424b87d60b3SJames Wright &sgs_dd_data->op_nodal_dd_inputs_ctx)); 425b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs)); 426b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs)); 427b87d60b3SJames Wright } 428b87d60b3SJames Wright 429b87d60b3SJames Wright { // Create operator for data-driven output handling 430b87d60b3SJames Wright CeedQFunction qf_sgs_dd_outputs; 431b87d60b3SJames Wright CeedOperator op_sgs_dd_outputs; 432b87d60b3SJames Wright 433b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc, 434b87d60b3SJames Wright &qf_sgs_dd_outputs)); 435b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx)); 436b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 437b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 438b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 439b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 440b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 441b87d60b3SJames Wright 442b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs)); 443b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 444b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 445b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 446b87d60b3SJames Wright PetscCallCeed(ceed, 447b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 448b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 449b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 450b87d60b3SJames Wright 451b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_outputs, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_outputs, NULL, sgs_dd_data->sgs_nodal_ceed, 452b87d60b3SJames Wright NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx)); 453b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs)); 454b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs)); 455b87d60b3SJames Wright } 456b87d60b3SJames Wright 457b87d60b3SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential; 4584c07ec22SJames Wright 4594c07ec22SJames Wright if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) { 4604c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed; 4614c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 462b87d60b3SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4634c07ec22SJames Wright } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) { 4644c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch; 4654c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 4664c07ec22SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4674c07ec22SJames Wright } 468b87d60b3SJames Wright 469b87d60b3SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 470b87d60b3SJames Wright 471b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 472b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&eigvec)); 473b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 474b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec)); 475b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs)); 476b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs)); 477fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 478b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 479b87d60b3SJames Wright } 480b87d60b3SJames Wright 4819c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual 482e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 48382baf964SJames Wright SgsDDData sgs_dd_data; 484be29160dSJames Wright CeedInt num_comp_q, q_data_size, num_comp_x; 485defe8520SJames Wright PetscInt dim; 4869c678832SJames Wright CeedQFunction qf_sgs_apply; 4879c678832SJames Wright CeedOperator op_sgs_apply; 4889c678832SJames Wright CeedBasis basis_sgs; 489be29160dSJames Wright CeedVector q_data; 490be29160dSJames Wright CeedElemRestriction elem_restr_qd; 4919c678832SJames Wright 4929c678832SJames Wright PetscFunctionBeginUser; 4930c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 4940c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 495e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 496e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 4979c678832SJames Wright 498be29160dSJames Wright { 499be29160dSJames Wright DMLabel domain_label = NULL; 500be29160dSJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 501be29160dSJames Wright 502be29160dSJames Wright PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs)); 503e3663b90SJames Wright PetscCall(QDataGet(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 504e3663b90SJames Wright &q_data, &q_data_size)); 505be29160dSJames Wright } 5069c678832SJames Wright 5070c373b74SJames Wright switch (honee->phys->state_var) { 5089c678832SJames Wright case STATEVAR_PRIMITIVE: 50942454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 5109c678832SJames Wright break; 5119c678832SJames Wright case STATEVAR_CONSERVATIVE: 51242454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 5139c678832SJames Wright break; 5149b103f75SJames Wright case STATEVAR_ENTROPY: 5159b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply)); 5169b103f75SJames Wright break; 5179c678832SJames Wright } 5189c678832SJames Wright 51940816385SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx)); 520b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 521be29160dSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE)); 522b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 523b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 5249c678832SJames Wright 525b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 526e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 527be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 528b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 529e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 5309c678832SJames Wright 5319c678832SJames Wright PetscCall( 5320c373b74SJames Wright OperatorApplyContextCreate(honee->dm, honee->dm, ceed, op_sgs_apply, honee->q_ceed, honee->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx)); 5339c678832SJames Wright 534be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 535be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 53687edb941SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs)); 537b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 538b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 539d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5409c678832SJames Wright } 5419c678832SJames Wright 5429c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual 5430c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) { 54482baf964SJames Wright SgsDDData sgs_dd_data; 5459c678832SJames Wright Vec VelocityGradient, SGSNodal_loc; 546cceb3143SJames Wright PetscMemType sgs_nodal_mem_type; 5478340219bSJames Wright NodalProjectionData grad_velo_proj; 5489c678832SJames Wright 5499c678832SJames Wright PetscFunctionBeginUser; 550ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 5510c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 5520c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 5538340219bSJames Wright PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 5548340219bSJames Wright PetscCall(VelocityGradientProjectionApply(grad_velo_proj, Q_loc, VelocityGradient)); 5559c678832SJames Wright 5569c678832SJames Wright // -- Compute Nodal SGS tensor 5579c678832SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 5580c373b74SJames Wright PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc)); 5599c678832SJames Wright 5609c678832SJames Wright // -- Compute contribution of the SGS stress 561a7dac1d5SJames Wright PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 5629c678832SJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 5639c678832SJames Wright 5649c678832SJames Wright // -- Return local SGS vector 565a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 5669c678832SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 5678340219bSJames Wright PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 568ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 569d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5709c678832SJames Wright } 5719c678832SJames Wright 57262b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN 573cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 57462b7942eSJames Wright PetscFunctionBeginUser; 57562b7942eSJames Wright for (PetscInt i = 0; i < N; i++) { 57662b7942eSJames Wright for (PetscInt j = 0; j < M; j++) { 57762b7942eSJames Wright B[j * N + i] = A[i * M + j]; 57862b7942eSJames Wright } 57962b7942eSJames Wright } 580d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 58162b7942eSJames Wright } 58262b7942eSJames Wright 58362b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct 584ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) { 585ad494f68SJames Wright SgsDDContext sgsdd_ctx; 58662b7942eSJames Wright PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 58762b7942eSJames Wright char file_path[PETSC_MAX_PATH_LEN]; 58862b7942eSJames Wright PetscScalar *temp; 58962b7942eSJames Wright 59062b7942eSJames Wright PetscFunctionBeginUser; 59162b7942eSJames Wright { 592ad494f68SJames Wright SgsDDContext sgsdd_temp; 59362b7942eSJames Wright PetscCall(PetscNew(&sgsdd_temp)); 59462b7942eSJames Wright *sgsdd_temp = **psgsdd_ctx; 59562b7942eSJames Wright sgsdd_temp->offsets.bias1 = 0; 59662b7942eSJames Wright sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 59762b7942eSJames Wright sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 59862b7942eSJames Wright sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 59962b7942eSJames Wright sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 60062b7942eSJames Wright PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 60162b7942eSJames Wright sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 60262b7942eSJames Wright PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 60362b7942eSJames Wright *sgsdd_ctx = *sgsdd_temp; 60462b7942eSJames Wright PetscCall(PetscFree(sgsdd_temp)); 60562b7942eSJames Wright } 60662b7942eSJames Wright 60762b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 60842454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 60962b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 61042454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 61162b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 61242454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 61362b7942eSJames Wright 61462b7942eSJames Wright { 61562b7942eSJames Wright PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 61662b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 61742454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 61862b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 61962b7942eSJames Wright PetscCall(PetscFree(temp)); 62062b7942eSJames Wright } 62162b7942eSJames Wright { 62262b7942eSJames Wright PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 62362b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 62442454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 62562b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 62662b7942eSJames Wright PetscCall(PetscFree(temp)); 62762b7942eSJames Wright } 62862b7942eSJames Wright 62962b7942eSJames Wright PetscCall(PetscFree(*psgsdd_ctx)); 63062b7942eSJames Wright *psgsdd_ctx = sgsdd_ctx; 631d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63262b7942eSJames Wright } 63362b7942eSJames Wright 634e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) { 635ee1455b7SJames Wright PetscReal alpha = 0; 636ad494f68SJames Wright SgsDDContext sgsdd_ctx; 6370c373b74SJames Wright MPI_Comm comm = honee->comm; 638ee1455b7SJames Wright char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 639ad494f68SJames Wright SgsDDSetupData sgs_dd_setup_data; 640cde3d787SJames Wright NewtonianIdealGasContext newt_ctx; 6418340219bSJames Wright NodalProjectionData grad_velo_proj; 64282baf964SJames Wright SgsDDData sgs_dd_data; 64362b7942eSJames Wright 64406f41313SJames Wright PetscFunctionBeginUser; 6458340219bSJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &grad_velo_proj)); 6460c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, GRAD_VELO_PROJ_KEY, grad_velo_proj, (PetscCtxDestroyFn *)NodalProjectionDataDestroy)); 6479ab09d51SJames Wright 64882baf964SJames Wright PetscCall(PetscNew(&sgs_dd_data)); 64982baf964SJames Wright sgs_dd_data->num_comp_inputs = 6; 65082baf964SJames Wright sgs_dd_data->num_comp_outputs = 6; 6510c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, SGS_DD_DATA_KEY, sgs_dd_data, (PetscCtxDestroyFn *)SgsDDDataDestroy)); 65262b7942eSJames Wright 6534c07ec22SJames Wright PetscCall(PetscNew(&sgs_dd_setup_data)); 6544c07ec22SJames Wright 6554b0f6111SJames Wright PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 65662b7942eSJames Wright PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 65762b7942eSJames Wright PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 65862b7942eSJames Wright sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 6594c07ec22SJames Wright PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead")); 6604c07ec22SJames Wright sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED; 6614c07ec22SJames Wright PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations, 6624c07ec22SJames Wright (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation, 6634c07ec22SJames Wright NULL)); 66462b7942eSJames Wright PetscOptionsEnd(); 66562b7942eSJames Wright 666b87d60b3SJames Wright PetscCall(PetscNew(&sgsdd_ctx)); 667f5dc303cSJames Wright *sgsdd_ctx = (struct SgsDDContext_){ 668f5dc303cSJames Wright .num_layers = 1, 669f5dc303cSJames Wright .num_inputs = 6, 670f5dc303cSJames Wright .num_outputs = 6, 671f5dc303cSJames Wright .num_neurons = 20, 672f5dc303cSJames Wright .alpha = alpha, 673f5dc303cSJames Wright }; 67462b7942eSJames Wright 675ad494f68SJames Wright PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 67662b7942eSJames Wright 677ee1455b7SJames Wright // -- Create DM for storing SGS tensor at nodes 67882baf964SJames Wright PetscCall(SgsDDCreateDM(honee->dm, &sgs_dd_data->dm_sgs, honee->app_ctx->degree, honee->app_ctx->q_extra, &sgs_dd_data->num_comp_sgs)); 679ee1455b7SJames Wright 680cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx)); 681cde3d787SJames Wright sgsdd_ctx->newt_ctx = *newt_ctx; 682cde3d787SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx)); 6830c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 684b4c37c5cSJames Wright PetscCallCeed(ceed, 685b4c37c5cSJames Wright CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 686b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 687ee1455b7SJames Wright 688e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx)); 68940816385SJames Wright 690c38c977aSJames Wright // -- Compute and store anisotropy tensor 691e3663b90SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed)); 692c38c977aSJames Wright 693ee1455b7SJames Wright // -- Create Nodal Evaluation Operator 6944c07ec22SJames Wright switch (sgs_dd_setup_data->sgs_dd_model_implementation) { 6954c07ec22SJames Wright case SGS_MODEL_DD_FUSED: 696e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data)); 6974c07ec22SJames Wright break; 6984c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_CEED: 6994c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_TORCH: 700e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data)); 7014c07ec22SJames Wright break; 7024c07ec22SJames Wright } 703ee1455b7SJames Wright 7049c678832SJames Wright // -- Create Operator to evalutate residual of SGS stress 705e3663b90SJames Wright PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data)); 7069c678832SJames Wright 707ad494f68SJames Wright PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data)); 708d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 70962b7942eSJames Wright } 710