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; 2782baf964SJames Wright PetscFunctionBeginUser; 2882baf964SJames Wright if (!sgs_dd_data_) PetscFunctionReturn(PETSC_SUCCESS); 2982baf964SJames Wright Ceed ceed = sgs_dd_data_->op_sgs_apply_ctx->ceed; 3082baf964SJames Wright 3182baf964SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->sgs_nodal_ceed)); 3282baf964SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->grad_velo_ceed)); 3382baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_evaluation_ctx)); 3482baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_sgs_apply_ctx)); 3582baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_inputs_ctx)); 3682baf964SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_outputs_ctx)); 3782baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_sgs)); 3882baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_inputs)); 3982baf964SJames Wright PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_outputs)); 4082baf964SJames 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)); 4182baf964SJames Wright PetscCall(PetscFree(sgs_dd_data_)); 4282baf964SJames Wright *sgs_dd_data = NULL; 4382baf964SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4482baf964SJames Wright } 4582baf964SJames Wright 46c38c977aSJames Wright typedef struct { 47c38c977aSJames Wright CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 48c38c977aSJames Wright CeedVector grid_aniso_ceed; 4940816385SJames Wright CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx; 504c07ec22SJames Wright SGSModelDDImplementation sgs_dd_model_implementation; 51ad494f68SJames Wright } *SgsDDSetupData; 52c38c977aSJames Wright 538340219bSJames Wright #define GRAD_VELO_PROJ_KEY "Gradient of Velocity Projection" 5482baf964SJames Wright #define SGS_DD_DATA_KEY "SGS Data Driven Data" 558340219bSJames Wright 56ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) { 57b4c37c5cSJames Wright Ceed ceed; 58ad494f68SJames Wright 59c38c977aSJames Wright PetscFunctionBeginUser; 60b4c37c5cSJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 61ad494f68SJames Wright 62b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 63b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 64b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 65b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 6640816385SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx)); 67c38c977aSJames Wright PetscCall(PetscFree(sgs_dd_setup_data)); 68519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 69d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 70c38c977aSJames Wright } 71c38c977aSJames Wright 72ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes 73ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 74ee1455b7SJames Wright PetscSection section; 75ee1455b7SJames Wright 76ee1455b7SJames Wright PetscFunctionBeginUser; 77ee1455b7SJames Wright *num_components = 6; 78ee1455b7SJames Wright 79ee1455b7SJames Wright PetscCall(DMClone(dm_source, dm_sgs)); 800dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE)); 81ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 82ee1455b7SJames Wright 83da4ca0cfSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 84ee1455b7SJames Wright 85ee1455b7SJames Wright PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 86ee1455b7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 87ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 88ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 89ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 90ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 91ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 92ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 93d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 94ee1455b7SJames Wright }; 95ee1455b7SJames Wright 96ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method 970c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 9882baf964SJames Wright SgsDDData sgs_dd_data; 99cceb3143SJames Wright PetscMemType q_mem_type; 100cceb3143SJames Wright 101cceb3143SJames Wright PetscFunctionBeginUser; 102*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 1030c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 104cceb3143SJames Wright 105cceb3143SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 106cceb3143SJames Wright 1070c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 108cceb3143SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 109cceb3143SJames Wright } 110cceb3143SJames Wright 111b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator 112e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 11382baf964SJames Wright SgsDDData sgs_dd_data; 1145930f037SJames Wright CeedQFunction qf_sgs_dd_nodal; 1155930f037SJames Wright CeedOperator op_sgs_dd_nodal; 11615c18037SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 117defe8520SJames Wright PetscInt dim; 1185930f037SJames Wright CeedVector inv_multiplicity; 119ee1455b7SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 12015c18037SJames Wright DMLabel domain_label = NULL; 12115c18037SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 1228340219bSJames Wright NodalProjectionData grad_velo_proj; 123ee1455b7SJames Wright 124ee1455b7SJames Wright PetscFunctionBeginUser; 1250c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 126e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 127e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 128b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 129*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 130*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 131ee1455b7SJames Wright 132ee1455b7SJames Wright { // Get velocity gradient information 133ee1455b7SJames Wright CeedOperatorField op_field; 1348340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 135b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 137ee1455b7SJames Wright } 13815c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 139b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 140ee1455b7SJames Wright 1415930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 1425930f037SJames Wright &inv_multiplicity)); 143ee1455b7SJames Wright 144ee1455b7SJames Wright // -- Create operator for SGS DD model nodal evaluation 1450c373b74SJames Wright switch (honee->phys->state_var) { 146ee1455b7SJames Wright case STATEVAR_PRIMITIVE: 147ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal)); 148ee1455b7SJames Wright break; 149ee1455b7SJames Wright case STATEVAR_CONSERVATIVE: 150ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal)); 151ee1455b7SJames Wright break; 1529b103f75SJames Wright case STATEVAR_ENTROPY: 1539b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal)); 1549b103f75SJames Wright break; 155ee1455b7SJames Wright } 156ee1455b7SJames Wright 157ee1455b7SJames Wright // Mesh/geometry order and solution basis order may differ, therefore must interpolate 158ee1455b7SJames Wright CeedBasis basis_x_to_q; 159e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q)); 160ee1455b7SJames Wright 161b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 162b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 163b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 164b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 165b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 166b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 167b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 168ee1455b7SJames Wright 169b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 170e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 171e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord)); 17258e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 17358e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 174b4c37c5cSJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 17558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 17658e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 177ee1455b7SJames Wright 1788340219bSJames 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, 1798340219bSJames Wright &sgs_dd_data->op_nodal_evaluation_ctx)); 180ee1455b7SJames Wright 181ee1455b7SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 182ad494f68SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Fused; 183ee1455b7SJames Wright 184b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 185b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 186b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 187fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 188b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 189b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 190d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 191ee1455b7SJames Wright } 192ee1455b7SJames Wright 1934c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation 1944c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 1954c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 196b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 197b87d60b3SJames Wright void **ctx) { 198b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inference; 199b87d60b3SJames Wright CeedOperator op_sgs_dd_inference; 200b87d60b3SJames Wright OperatorApplyContext *op_context = (OperatorApplyContext *)ctx; 201b87d60b3SJames Wright 202b87d60b3SJames Wright PetscFunctionBeginUser; 203b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc, 204b87d60b3SJames Wright &qf_sgs_dd_inference)); 205b87d60b3SJames Wright 206b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx)); 207b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 208b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE)); 209b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 210b87d60b3SJames Wright 211b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference)); 212b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 213b87d60b3SJames Wright PetscCallCeed(ceed, 214b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 215b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 216b87d60b3SJames Wright 217b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL, 218b87d60b3SJames Wright op_context)); 219b87d60b3SJames Wright sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy; 220b87d60b3SJames Wright 221b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference)); 222b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference)); 223b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 224b87d60b3SJames Wright } 225b87d60b3SJames Wright 2264c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation 2274c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 228b87d60b3SJames Wright OperatorApplyContext op_context = *(OperatorApplyContext *)ctx; 229b87d60b3SJames Wright 230b87d60b3SJames Wright PetscFunctionBeginUser; 231ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 232ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 233b40a7e63SJames Wright PetscCall(PetscLogGpuTimeBegin()); 234b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context)); 235b40a7e63SJames Wright PetscCall(PetscLogGpuTimeEnd()); 236ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 237ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 238b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 239b87d60b3SJames Wright } 240b87d60b3SJames Wright 2414c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch 2424c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 2434c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 2444c07ec22SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 2454c07ec22SJames Wright void **ctx) { 2464c07ec22SJames Wright const char *ceed_resource; 2474c07ec22SJames Wright char model_path[PETSC_MAX_PATH_LEN] = ""; 2484c07ec22SJames Wright TorchDeviceType model_device_type; 2494c07ec22SJames Wright 2504c07ec22SJames Wright PetscFunctionBeginUser; 2514c07ec22SJames Wright PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource)); 2524c07ec22SJames Wright if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA; 2534c07ec22SJames Wright else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP; 2547ffa0ff8SJames Wright // On-device XPU is not working reliably currently, default to CPU inference evaluation 2557ffa0ff8SJames Wright // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU; 2564c07ec22SJames Wright else model_device_type = TORCH_DEVICE_CPU; 2574c07ec22SJames Wright PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL)); 2584c07ec22SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL)); 2594c07ec22SJames Wright 2604c07ec22SJames Wright PetscCall(LoadModel_Torch(model_path, model_device_type)); 2614c07ec22SJames Wright 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; 2684c07ec22SJames Wright PetscFunctionBeginUser; 2694c07ec22SJames Wright if (!run_through) { 2704c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view")); 2714c07ec22SJames Wright } 2724c07ec22SJames Wright PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc)); 2734c07ec22SJames Wright if (!run_through) { 2744c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view")); 2754c07ec22SJames Wright run_through = PETSC_TRUE; 2764c07ec22SJames Wright } 2774c07ec22SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2784c07ec22SJames Wright } 2794c07ec22SJames Wright 280b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method 2810c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 28282baf964SJames Wright SgsDDData sgs_dd_data; 283b87d60b3SJames Wright PetscMemType q_mem_type; 284b87d60b3SJames Wright Vec DD_Inputs_loc, DD_Outputs_loc; 285b87d60b3SJames Wright 286b87d60b3SJames Wright PetscFunctionBeginUser; 287*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 288b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 289b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 2900c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 291b87d60b3SJames Wright 292b87d60b3SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx)); 293b87d60b3SJames Wright PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx)); 294b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx)); 295b87d60b3SJames Wright 2960c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 297b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 298b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 299b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 300b87d60b3SJames Wright } 301b87d60b3SJames Wright 302b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators 303e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 30482baf964SJames Wright SgsDDData sgs_dd_data; 305b87d60b3SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1; 306b87d60b3SJames Wright PetscInt dim; 307b87d60b3SJames Wright CeedVector inv_multiplicity, eigvec; 3088340219bSJames Wright NodalProjectionData grad_velo_proj; 309b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs, 310b87d60b3SJames Wright elem_restr_dd_outputs; 311b87d60b3SJames Wright DMLabel domain_label = NULL; 312b87d60b3SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 313b87d60b3SJames Wright 314b87d60b3SJames Wright PetscFunctionBeginUser; 315*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 316b87d60b3SJames Wright { // Create DMs for data-driven input and output values 317b87d60b3SJames Wright PetscSection section; 318b87d60b3SJames Wright PetscInt degree, q_extra; 319b87d60b3SJames Wright { // Get degree and number of quadrature points from dm_sgs 320b87d60b3SJames Wright PetscFE fe; 321b87d60b3SJames Wright PetscSpace basis; 322b87d60b3SJames Wright PetscQuadrature quadrature; 323b87d60b3SJames Wright PetscInt num_qpnts; 324b87d60b3SJames Wright PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe)); 325b87d60b3SJames Wright PetscCall(PetscFEGetBasisSpace(fe, &basis)); 326b87d60b3SJames Wright PetscCall(PetscSpaceGetDegree(basis, °ree, NULL)); 327b87d60b3SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 328b87d60b3SJames Wright PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts)); 329b87d60b3SJames Wright q_extra = degree - num_qpnts; 330b87d60b3SJames Wright } 331b87d60b3SJames Wright 332b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs)); 3330dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE)); 334b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs")); 335b87d60b3SJames 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)); 336b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, §ion)); 337b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 338b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) { 339b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 340b87d60b3SJames Wright 341b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1)); 342b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 343b87d60b3SJames Wright } 344b87d60b3SJames Wright 345b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs)); 3460dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE)); 347b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs")); 348b87d60b3SJames 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)); 349b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, §ion)); 350b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 351b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) { 352b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 353b87d60b3SJames Wright 354b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1)); 355b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 356b87d60b3SJames Wright } 357b87d60b3SJames Wright } 358b87d60b3SJames Wright 3590c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 360e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 361e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 362b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 363*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 364b87d60b3SJames Wright 365b87d60b3SJames Wright { // Get velocity gradient information 366b87d60b3SJames Wright CeedOperatorField op_field; 3678340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 368b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 369b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 370b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL)); 371b87d60b3SJames Wright } 372b87d60b3SJames Wright 373b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 374b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 375b87d60b3SJames Wright PetscCall( 376b87d60b3SJames Wright DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec)); 377b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL)); 378b87d60b3SJames Wright 379b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs)); 380b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs)); 381b87d60b3SJames Wright 3825930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 3835930f037SJames Wright &inv_multiplicity)); 384b87d60b3SJames Wright 385b87d60b3SJames Wright { // Create operator for data-driven input evaluation 386b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inputs; 387b87d60b3SJames Wright CeedOperator op_sgs_dd_inputs; 388b87d60b3SJames Wright 3890c373b74SJames Wright switch (honee->phys->state_var) { 390b87d60b3SJames Wright case STATEVAR_PRIMITIVE: 391b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim, 392b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs)); 393b87d60b3SJames Wright break; 394b87d60b3SJames Wright case STATEVAR_CONSERVATIVE: 395b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv, 396b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs)); 397b87d60b3SJames Wright break; 3989b103f75SJames Wright case STATEVAR_ENTROPY: 3999b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy, 4009b103f75SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs)); 4019b103f75SJames Wright break; 402b87d60b3SJames Wright } 403b87d60b3SJames Wright 404b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx)); 405b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE)); 406b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 407b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 408b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 409b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 410b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 411b87d60b3SJames Wright 412b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs)); 413e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 414b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 415b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 416b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 417b87d60b3SJames Wright PetscCallCeed(ceed, 418b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 419b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 420b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 421b87d60b3SJames Wright 4228340219bSJames Wright PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL, 423b87d60b3SJames Wright &sgs_dd_data->op_nodal_dd_inputs_ctx)); 424b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs)); 425b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs)); 426b87d60b3SJames Wright } 427b87d60b3SJames Wright 428b87d60b3SJames Wright { // Create operator for data-driven output handling 429b87d60b3SJames Wright CeedQFunction qf_sgs_dd_outputs; 430b87d60b3SJames Wright CeedOperator op_sgs_dd_outputs; 431b87d60b3SJames Wright 432b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc, 433b87d60b3SJames Wright &qf_sgs_dd_outputs)); 434b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx)); 435b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 436b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 437b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 438b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 439b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 440b87d60b3SJames Wright 441b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs)); 442b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 443b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 444b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 445b87d60b3SJames Wright PetscCallCeed(ceed, 446b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 447b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 448b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 449b87d60b3SJames Wright 450b87d60b3SJames 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, 451b87d60b3SJames Wright NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx)); 452b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs)); 453b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs)); 454b87d60b3SJames Wright } 455b87d60b3SJames Wright 456b87d60b3SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential; 4574c07ec22SJames Wright 4584c07ec22SJames Wright if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) { 4594c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed; 4604c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 461b87d60b3SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4624c07ec22SJames Wright } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) { 4634c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch; 4644c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 4654c07ec22SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4664c07ec22SJames Wright } 467b87d60b3SJames Wright 468b87d60b3SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 469b87d60b3SJames Wright 470b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 471b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&eigvec)); 472b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 473b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec)); 474b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs)); 475b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs)); 476fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 477b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 478b87d60b3SJames Wright } 479b87d60b3SJames Wright 4809c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual 481e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 48282baf964SJames Wright SgsDDData sgs_dd_data; 483be29160dSJames Wright CeedInt num_comp_q, q_data_size, num_comp_x; 484defe8520SJames Wright PetscInt dim; 4859c678832SJames Wright CeedQFunction qf_sgs_apply; 4869c678832SJames Wright CeedOperator op_sgs_apply; 4879c678832SJames Wright CeedBasis basis_sgs; 488be29160dSJames Wright CeedVector q_data; 489be29160dSJames Wright CeedElemRestriction elem_restr_qd; 4909c678832SJames Wright 4919c678832SJames Wright PetscFunctionBeginUser; 492*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 4930c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 494e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 495e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 4969c678832SJames Wright 497be29160dSJames Wright { 498be29160dSJames Wright DMLabel domain_label = NULL; 499be29160dSJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 500be29160dSJames Wright 501be29160dSJames Wright PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs)); 502e3663b90SJames 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, 503e3663b90SJames Wright &q_data, &q_data_size)); 504be29160dSJames Wright } 5059c678832SJames Wright 5060c373b74SJames Wright switch (honee->phys->state_var) { 5079c678832SJames Wright case STATEVAR_PRIMITIVE: 50842454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 5099c678832SJames Wright break; 5109c678832SJames Wright case STATEVAR_CONSERVATIVE: 51142454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 5129c678832SJames Wright break; 5139b103f75SJames Wright case STATEVAR_ENTROPY: 5149b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply)); 5159b103f75SJames Wright break; 5169c678832SJames Wright } 5179c678832SJames Wright 51840816385SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx)); 519b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 520be29160dSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE)); 521b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 522b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 5239c678832SJames Wright 524b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 525e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 526be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 527b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 528e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 5299c678832SJames Wright 5309c678832SJames Wright PetscCall( 5310c373b74SJames 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)); 5329c678832SJames Wright 533be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 534be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 53587edb941SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs)); 536b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 537b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 538d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5399c678832SJames Wright } 5409c678832SJames Wright 5419c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual 5420c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) { 54382baf964SJames Wright SgsDDData sgs_dd_data; 5449c678832SJames Wright Vec VelocityGradient, SGSNodal_loc; 545cceb3143SJames Wright PetscMemType sgs_nodal_mem_type; 5468340219bSJames Wright NodalProjectionData grad_velo_proj; 5479c678832SJames Wright 5489c678832SJames Wright PetscFunctionBeginUser; 549ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 550*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 551*0c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 5528340219bSJames Wright PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 5538340219bSJames Wright PetscCall(VelocityGradientProjectionApply(grad_velo_proj, Q_loc, VelocityGradient)); 5549c678832SJames Wright 5559c678832SJames Wright // -- Compute Nodal SGS tensor 5569c678832SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 5570c373b74SJames Wright PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc)); 5589c678832SJames Wright 5599c678832SJames Wright // -- Compute contribution of the SGS stress 560a7dac1d5SJames Wright PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 5619c678832SJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 5629c678832SJames Wright 5639c678832SJames Wright // -- Return local SGS vector 564a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 5659c678832SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 5668340219bSJames Wright PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 567ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 568d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5699c678832SJames Wright } 5709c678832SJames Wright 57162b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN 572cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 57362b7942eSJames Wright PetscFunctionBeginUser; 57462b7942eSJames Wright for (PetscInt i = 0; i < N; i++) { 57562b7942eSJames Wright for (PetscInt j = 0; j < M; j++) { 57662b7942eSJames Wright B[j * N + i] = A[i * M + j]; 57762b7942eSJames Wright } 57862b7942eSJames Wright } 579d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 58062b7942eSJames Wright } 58162b7942eSJames Wright 58262b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct 583ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) { 584ad494f68SJames Wright SgsDDContext sgsdd_ctx; 58562b7942eSJames Wright PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 58662b7942eSJames Wright char file_path[PETSC_MAX_PATH_LEN]; 58762b7942eSJames Wright PetscScalar *temp; 58862b7942eSJames Wright 58962b7942eSJames Wright PetscFunctionBeginUser; 59062b7942eSJames Wright { 591ad494f68SJames Wright SgsDDContext sgsdd_temp; 59262b7942eSJames Wright PetscCall(PetscNew(&sgsdd_temp)); 59362b7942eSJames Wright *sgsdd_temp = **psgsdd_ctx; 59462b7942eSJames Wright sgsdd_temp->offsets.bias1 = 0; 59562b7942eSJames Wright sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 59662b7942eSJames Wright sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 59762b7942eSJames Wright sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 59862b7942eSJames Wright sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 59962b7942eSJames Wright PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 60062b7942eSJames Wright sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 60162b7942eSJames Wright PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 60262b7942eSJames Wright *sgsdd_ctx = *sgsdd_temp; 60362b7942eSJames Wright PetscCall(PetscFree(sgsdd_temp)); 60462b7942eSJames Wright } 60562b7942eSJames Wright 60662b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 60742454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 60862b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 60942454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 61062b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 61142454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 61262b7942eSJames Wright 61362b7942eSJames Wright { 61462b7942eSJames Wright PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 61562b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 61642454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 61762b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 61862b7942eSJames Wright PetscCall(PetscFree(temp)); 61962b7942eSJames Wright } 62062b7942eSJames Wright { 62162b7942eSJames Wright PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 62262b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 62342454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 62462b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 62562b7942eSJames Wright PetscCall(PetscFree(temp)); 62662b7942eSJames Wright } 62762b7942eSJames Wright 62862b7942eSJames Wright PetscCall(PetscFree(*psgsdd_ctx)); 62962b7942eSJames Wright *psgsdd_ctx = sgsdd_ctx; 630d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63162b7942eSJames Wright } 63262b7942eSJames Wright 633e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) { 634ee1455b7SJames Wright PetscReal alpha = 0; 635ad494f68SJames Wright SgsDDContext sgsdd_ctx; 6360c373b74SJames Wright MPI_Comm comm = honee->comm; 637ee1455b7SJames Wright char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 638ad494f68SJames Wright SgsDDSetupData sgs_dd_setup_data; 639ee1455b7SJames Wright NewtonianIdealGasContext gas; 6408340219bSJames Wright NodalProjectionData grad_velo_proj; 64182baf964SJames Wright SgsDDData sgs_dd_data; 64262b7942eSJames Wright 64306f41313SJames Wright PetscFunctionBeginUser; 6448340219bSJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &grad_velo_proj)); 645*0c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, GRAD_VELO_PROJ_KEY, grad_velo_proj, (PetscCtxDestroyFn *)NodalProjectionDataDestroy)); 6469ab09d51SJames Wright 64782baf964SJames Wright PetscCall(PetscNew(&sgs_dd_data)); 64882baf964SJames Wright sgs_dd_data->num_comp_inputs = 6; 64982baf964SJames Wright sgs_dd_data->num_comp_outputs = 6; 650*0c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, SGS_DD_DATA_KEY, sgs_dd_data, (PetscCtxDestroyFn *)SgsDDDataDestroy)); 65162b7942eSJames Wright 6524c07ec22SJames Wright PetscCall(PetscNew(&sgs_dd_setup_data)); 6534c07ec22SJames Wright 6544b0f6111SJames Wright PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 65562b7942eSJames Wright PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 65662b7942eSJames Wright PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 65762b7942eSJames Wright sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 6584c07ec22SJames Wright PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead")); 6594c07ec22SJames Wright sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED; 6604c07ec22SJames Wright PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations, 6614c07ec22SJames Wright (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation, 6624c07ec22SJames Wright NULL)); 66362b7942eSJames Wright PetscOptionsEnd(); 66462b7942eSJames Wright 665b87d60b3SJames Wright PetscCall(PetscNew(&sgsdd_ctx)); 6664b0f6111SJames Wright sgsdd_ctx->num_layers = 1; 66762b7942eSJames Wright sgsdd_ctx->num_inputs = 6; 66862b7942eSJames Wright sgsdd_ctx->num_outputs = 6; 66962b7942eSJames Wright sgsdd_ctx->num_neurons = 20; 67062b7942eSJames Wright sgsdd_ctx->alpha = alpha; 67162b7942eSJames Wright 672ad494f68SJames Wright PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 67362b7942eSJames Wright 674ee1455b7SJames Wright // -- Create DM for storing SGS tensor at nodes 67582baf964SJames 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)); 676ee1455b7SJames Wright 677e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas)); 678ee1455b7SJames Wright sgsdd_ctx->gas = *gas; 679e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 6800c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 681b4c37c5cSJames Wright PetscCallCeed(ceed, 682b4c37c5cSJames Wright CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 683b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 684ee1455b7SJames Wright 685e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx)); 68640816385SJames Wright 687c38c977aSJames Wright // -- Compute and store anisotropy tensor 688e3663b90SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed)); 689c38c977aSJames Wright 690ee1455b7SJames Wright // -- Create Nodal Evaluation Operator 6914c07ec22SJames Wright switch (sgs_dd_setup_data->sgs_dd_model_implementation) { 6924c07ec22SJames Wright case SGS_MODEL_DD_FUSED: 693e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data)); 6944c07ec22SJames Wright break; 6954c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_CEED: 6964c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_TORCH: 697e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data)); 6984c07ec22SJames Wright break; 6994c07ec22SJames Wright } 700ee1455b7SJames Wright 7019c678832SJames Wright // -- Create Operator to evalutate residual of SGS stress 702e3663b90SJames Wright PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data)); 7039c678832SJames Wright 704ad494f68SJames Wright PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data)); 705d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 70662b7942eSJames Wright } 707