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 11c38c977aSJames Wright typedef struct { 12c38c977aSJames Wright CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 13c38c977aSJames Wright CeedVector grid_aniso_ceed; 1440816385SJames Wright CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx; 154c07ec22SJames Wright SGSModelDDImplementation sgs_dd_model_implementation; 16ad494f68SJames Wright } *SgsDDSetupData; 17c38c977aSJames Wright 18*8340219bSJames Wright #define GRAD_VELO_PROJ_KEY "Gradient of Velocity Projection" 19*8340219bSJames Wright 20ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) { 21b4c37c5cSJames Wright Ceed ceed; 22ad494f68SJames Wright 23c38c977aSJames Wright PetscFunctionBeginUser; 24b4c37c5cSJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 25ad494f68SJames Wright 26b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 27b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 28b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 29b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 3040816385SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx)); 31c38c977aSJames Wright PetscCall(PetscFree(sgs_dd_setup_data)); 32519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 33d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 34c38c977aSJames Wright } 35c38c977aSJames Wright 36ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes 37ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 38ee1455b7SJames Wright PetscSection section; 39ee1455b7SJames Wright 40ee1455b7SJames Wright PetscFunctionBeginUser; 41ee1455b7SJames Wright *num_components = 6; 42ee1455b7SJames Wright 43ee1455b7SJames Wright PetscCall(DMClone(dm_source, dm_sgs)); 440dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE)); 45ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 46ee1455b7SJames Wright 47da4ca0cfSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 48ee1455b7SJames Wright 49ee1455b7SJames Wright PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 50ee1455b7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 51ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 52ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 53ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 54ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 55ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 56ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 57d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 58ee1455b7SJames Wright }; 59ee1455b7SJames Wright 60ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method 610c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 620c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 63cceb3143SJames Wright PetscMemType q_mem_type; 64cceb3143SJames Wright 65cceb3143SJames Wright PetscFunctionBeginUser; 660c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 67cceb3143SJames Wright 68cceb3143SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 69cceb3143SJames Wright 700c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 71cceb3143SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 72cceb3143SJames Wright } 73cceb3143SJames Wright 74b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator 75e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 760c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 775930f037SJames Wright CeedQFunction qf_sgs_dd_nodal; 785930f037SJames Wright CeedOperator op_sgs_dd_nodal; 7915c18037SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 80defe8520SJames Wright PetscInt dim; 815930f037SJames Wright CeedVector inv_multiplicity; 82ee1455b7SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 8315c18037SJames Wright DMLabel domain_label = NULL; 8415c18037SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 85*8340219bSJames Wright NodalProjectionData grad_velo_proj; 86ee1455b7SJames Wright 87ee1455b7SJames Wright PetscFunctionBeginUser; 880c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 89e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 90e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 91b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 92*8340219bSJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 93ee1455b7SJames Wright 94ee1455b7SJames Wright { // Get velocity gradient information 95ee1455b7SJames Wright CeedOperatorField op_field; 96*8340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 97b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 98b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 99ee1455b7SJames Wright } 10015c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 101b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 102ee1455b7SJames Wright 1035930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 1045930f037SJames Wright &inv_multiplicity)); 105ee1455b7SJames Wright 106ee1455b7SJames Wright // -- Create operator for SGS DD model nodal evaluation 1070c373b74SJames Wright switch (honee->phys->state_var) { 108ee1455b7SJames Wright case STATEVAR_PRIMITIVE: 109ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal)); 110ee1455b7SJames Wright break; 111ee1455b7SJames Wright case STATEVAR_CONSERVATIVE: 112ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal)); 113ee1455b7SJames Wright break; 1149b103f75SJames Wright case STATEVAR_ENTROPY: 1159b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal)); 1169b103f75SJames Wright break; 117ee1455b7SJames Wright } 118ee1455b7SJames Wright 119ee1455b7SJames Wright // Mesh/geometry order and solution basis order may differ, therefore must interpolate 120ee1455b7SJames Wright CeedBasis basis_x_to_q; 121e3663b90SJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q)); 122ee1455b7SJames Wright 123b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 124b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 125b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 126b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 127b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 128b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 129b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 130ee1455b7SJames Wright 131b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 132e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 133e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord)); 13458e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 13558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 136b4c37c5cSJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 13758e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 13858e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 139ee1455b7SJames Wright 140*8340219bSJames 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, 141*8340219bSJames Wright &sgs_dd_data->op_nodal_evaluation_ctx)); 142ee1455b7SJames Wright 143ee1455b7SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 144ad494f68SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Fused; 145ee1455b7SJames Wright 146b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 147b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 148b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 149fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 150b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 151b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 152d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 153ee1455b7SJames Wright } 154ee1455b7SJames Wright 1554c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation 1564c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 1574c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 158b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 159b87d60b3SJames Wright void **ctx) { 160b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inference; 161b87d60b3SJames Wright CeedOperator op_sgs_dd_inference; 162b87d60b3SJames Wright OperatorApplyContext *op_context = (OperatorApplyContext *)ctx; 163b87d60b3SJames Wright 164b87d60b3SJames Wright PetscFunctionBeginUser; 165b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc, 166b87d60b3SJames Wright &qf_sgs_dd_inference)); 167b87d60b3SJames Wright 168b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx)); 169b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 170b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE)); 171b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 172b87d60b3SJames Wright 173b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference)); 174b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 175b87d60b3SJames Wright PetscCallCeed(ceed, 176b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 177b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 178b87d60b3SJames Wright 179b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL, 180b87d60b3SJames Wright op_context)); 181b87d60b3SJames Wright sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy; 182b87d60b3SJames Wright 183b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference)); 184b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference)); 185b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 186b87d60b3SJames Wright } 187b87d60b3SJames Wright 1884c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation 1894c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 190b87d60b3SJames Wright OperatorApplyContext op_context = *(OperatorApplyContext *)ctx; 191b87d60b3SJames Wright 192b87d60b3SJames Wright PetscFunctionBeginUser; 193ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 194ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 195b40a7e63SJames Wright PetscCall(PetscLogGpuTimeBegin()); 196b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context)); 197b40a7e63SJames Wright PetscCall(PetscLogGpuTimeEnd()); 198ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 199ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 200b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 201b87d60b3SJames Wright } 202b87d60b3SJames Wright 2034c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch 2044c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 2054c07ec22SJames Wright CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 2064c07ec22SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 2074c07ec22SJames Wright void **ctx) { 2084c07ec22SJames Wright const char *ceed_resource; 2094c07ec22SJames Wright char model_path[PETSC_MAX_PATH_LEN] = ""; 2104c07ec22SJames Wright TorchDeviceType model_device_type; 2114c07ec22SJames Wright 2124c07ec22SJames Wright PetscFunctionBeginUser; 2134c07ec22SJames Wright PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource)); 2144c07ec22SJames Wright if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA; 2154c07ec22SJames Wright else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP; 2167ffa0ff8SJames Wright // On-device XPU is not working reliably currently, default to CPU inference evaluation 2177ffa0ff8SJames Wright // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU; 2184c07ec22SJames Wright else model_device_type = TORCH_DEVICE_CPU; 2194c07ec22SJames Wright PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL)); 2204c07ec22SJames Wright PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL)); 2214c07ec22SJames Wright 2224c07ec22SJames Wright PetscCall(LoadModel_Torch(model_path, model_device_type)); 2234c07ec22SJames Wright 2244c07ec22SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2254c07ec22SJames Wright } 2264c07ec22SJames Wright 2274c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch 2284c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 2294c07ec22SJames Wright static PetscBool run_through = PETSC_FALSE; 2304c07ec22SJames Wright PetscFunctionBeginUser; 2314c07ec22SJames Wright if (!run_through) { 2324c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view")); 2334c07ec22SJames Wright } 2344c07ec22SJames Wright PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc)); 2354c07ec22SJames Wright if (!run_through) { 2364c07ec22SJames Wright PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view")); 2374c07ec22SJames Wright run_through = PETSC_TRUE; 2384c07ec22SJames Wright } 2394c07ec22SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 2404c07ec22SJames Wright } 2414c07ec22SJames Wright 242b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method 2430c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 2440c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 245b87d60b3SJames Wright PetscMemType q_mem_type; 246b87d60b3SJames Wright Vec DD_Inputs_loc, DD_Outputs_loc; 247b87d60b3SJames Wright 248b87d60b3SJames Wright PetscFunctionBeginUser; 249b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 250b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 2510c373b74SJames Wright PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 252b87d60b3SJames Wright 253b87d60b3SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx)); 254b87d60b3SJames Wright PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx)); 255b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx)); 256b87d60b3SJames Wright 2570c373b74SJames Wright PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 258b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 259b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 260b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 261b87d60b3SJames Wright } 262b87d60b3SJames Wright 263b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators 264e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 2650c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 266b87d60b3SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1; 267b87d60b3SJames Wright PetscInt dim; 268b87d60b3SJames Wright CeedVector inv_multiplicity, eigvec; 269*8340219bSJames Wright NodalProjectionData grad_velo_proj; 270b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs, 271b87d60b3SJames Wright elem_restr_dd_outputs; 272b87d60b3SJames Wright DMLabel domain_label = NULL; 273b87d60b3SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 274b87d60b3SJames Wright 275b87d60b3SJames Wright PetscFunctionBeginUser; 276b87d60b3SJames Wright { // Create DMs for data-driven input and output values 277b87d60b3SJames Wright PetscSection section; 278b87d60b3SJames Wright PetscInt degree, q_extra; 279b87d60b3SJames Wright { // Get degree and number of quadrature points from dm_sgs 280b87d60b3SJames Wright PetscFE fe; 281b87d60b3SJames Wright PetscSpace basis; 282b87d60b3SJames Wright PetscQuadrature quadrature; 283b87d60b3SJames Wright PetscInt num_qpnts; 284b87d60b3SJames Wright PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe)); 285b87d60b3SJames Wright PetscCall(PetscFEGetBasisSpace(fe, &basis)); 286b87d60b3SJames Wright PetscCall(PetscSpaceGetDegree(basis, °ree, NULL)); 287b87d60b3SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 288b87d60b3SJames Wright PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts)); 289b87d60b3SJames Wright q_extra = degree - num_qpnts; 290b87d60b3SJames Wright } 291b87d60b3SJames Wright 292b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs)); 2930dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE)); 294b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs")); 295b87d60b3SJames 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)); 296b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, §ion)); 297b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 298b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) { 299b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 300b87d60b3SJames Wright 301b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1)); 302b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 303b87d60b3SJames Wright } 304b87d60b3SJames Wright 305b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs)); 3060dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE)); 307b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs")); 308b87d60b3SJames 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)); 309b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, §ion)); 310b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 311b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) { 312b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 313b87d60b3SJames Wright 314b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1)); 315b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 316b87d60b3SJames Wright } 317b87d60b3SJames Wright } 318b87d60b3SJames Wright 3190c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 320e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 321e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 322b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 323*8340219bSJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 324b87d60b3SJames Wright 325b87d60b3SJames Wright { // Get velocity gradient information 326b87d60b3SJames Wright CeedOperatorField op_field; 327*8340219bSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 328b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 329b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 330b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL)); 331b87d60b3SJames Wright } 332b87d60b3SJames Wright 333b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 334b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 335b87d60b3SJames Wright PetscCall( 336b87d60b3SJames Wright DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec)); 337b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL)); 338b87d60b3SJames Wright 339b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs)); 340b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs)); 341b87d60b3SJames Wright 3425930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 3435930f037SJames Wright &inv_multiplicity)); 344b87d60b3SJames Wright 345b87d60b3SJames Wright { // Create operator for data-driven input evaluation 346b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inputs; 347b87d60b3SJames Wright CeedOperator op_sgs_dd_inputs; 348b87d60b3SJames Wright 3490c373b74SJames Wright switch (honee->phys->state_var) { 350b87d60b3SJames Wright case STATEVAR_PRIMITIVE: 351b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim, 352b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs)); 353b87d60b3SJames Wright break; 354b87d60b3SJames Wright case STATEVAR_CONSERVATIVE: 355b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv, 356b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs)); 357b87d60b3SJames Wright break; 3589b103f75SJames Wright case STATEVAR_ENTROPY: 3599b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy, 3609b103f75SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs)); 3619b103f75SJames Wright break; 362b87d60b3SJames Wright } 363b87d60b3SJames Wright 364b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx)); 365b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE)); 366b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 367b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 368b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 369b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 370b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 371b87d60b3SJames Wright 372b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs)); 373e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 374b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 375b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 376b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 377b87d60b3SJames Wright PetscCallCeed(ceed, 378b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 379b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 380b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 381b87d60b3SJames Wright 382*8340219bSJames Wright PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL, 383b87d60b3SJames Wright &sgs_dd_data->op_nodal_dd_inputs_ctx)); 384b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs)); 385b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs)); 386b87d60b3SJames Wright } 387b87d60b3SJames Wright 388b87d60b3SJames Wright { // Create operator for data-driven output handling 389b87d60b3SJames Wright CeedQFunction qf_sgs_dd_outputs; 390b87d60b3SJames Wright CeedOperator op_sgs_dd_outputs; 391b87d60b3SJames Wright 392b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc, 393b87d60b3SJames Wright &qf_sgs_dd_outputs)); 394b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx)); 395b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 396b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 397b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 398b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 399b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 400b87d60b3SJames Wright 401b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs)); 402b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 403b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 404b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 405b87d60b3SJames Wright PetscCallCeed(ceed, 406b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 407b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 408b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 409b87d60b3SJames Wright 410b87d60b3SJames 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, 411b87d60b3SJames Wright NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx)); 412b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs)); 413b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs)); 414b87d60b3SJames Wright } 415b87d60b3SJames Wright 416b87d60b3SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential; 4174c07ec22SJames Wright 4184c07ec22SJames Wright if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) { 4194c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed; 4204c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 421b87d60b3SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4224c07ec22SJames Wright } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) { 4234c07ec22SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch; 4244c07ec22SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 4254c07ec22SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 4264c07ec22SJames Wright } 427b87d60b3SJames Wright 428b87d60b3SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 429b87d60b3SJames Wright 430b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 431b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&eigvec)); 432b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 433b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec)); 434b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs)); 435b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs)); 436fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 437b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 438b87d60b3SJames Wright } 439b87d60b3SJames Wright 4409c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual 441e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 4420c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 443be29160dSJames Wright CeedInt num_comp_q, q_data_size, num_comp_x; 444defe8520SJames Wright PetscInt dim; 4459c678832SJames Wright CeedQFunction qf_sgs_apply; 4469c678832SJames Wright CeedOperator op_sgs_apply; 4479c678832SJames Wright CeedBasis basis_sgs; 448be29160dSJames Wright CeedVector q_data; 449be29160dSJames Wright CeedElemRestriction elem_restr_qd; 4509c678832SJames Wright 4519c678832SJames Wright PetscFunctionBeginUser; 4520c373b74SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 453e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 454e3663b90SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 4559c678832SJames Wright 456be29160dSJames Wright { 457be29160dSJames Wright DMLabel domain_label = NULL; 458be29160dSJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 459be29160dSJames Wright 460be29160dSJames Wright PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs)); 461e3663b90SJames 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, 462e3663b90SJames Wright &q_data, &q_data_size)); 463be29160dSJames Wright } 4649c678832SJames Wright 4650c373b74SJames Wright switch (honee->phys->state_var) { 4669c678832SJames Wright case STATEVAR_PRIMITIVE: 46742454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 4689c678832SJames Wright break; 4699c678832SJames Wright case STATEVAR_CONSERVATIVE: 47042454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 4719c678832SJames Wright break; 4729b103f75SJames Wright case STATEVAR_ENTROPY: 4739b103f75SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply)); 4749b103f75SJames Wright break; 4759c678832SJames Wright } 4769c678832SJames Wright 47740816385SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx)); 478b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 479be29160dSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE)); 480b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 481b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 4829c678832SJames Wright 483b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 484e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 485be29160dSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 486b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 487e3663b90SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 4889c678832SJames Wright 4899c678832SJames Wright PetscCall( 4900c373b74SJames 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)); 4919c678832SJames Wright 492be29160dSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 493be29160dSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 49487edb941SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs)); 495b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 496b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 497d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4989c678832SJames Wright } 4999c678832SJames Wright 5009c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual 5010c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) { 5020c373b74SJames Wright SgsDDData sgs_dd_data = honee->sgs_dd_data; 5039c678832SJames Wright Vec VelocityGradient, SGSNodal_loc; 504cceb3143SJames Wright PetscMemType sgs_nodal_mem_type; 505*8340219bSJames Wright NodalProjectionData grad_velo_proj; 5069c678832SJames Wright 5079c678832SJames Wright PetscFunctionBeginUser; 508ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 509*8340219bSJames Wright PetscCall(PetscObjectContainerQuery((PetscObject)honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 510*8340219bSJames Wright PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 511*8340219bSJames Wright PetscCall(VelocityGradientProjectionApply(grad_velo_proj, Q_loc, VelocityGradient)); 5129c678832SJames Wright 5139c678832SJames Wright // -- Compute Nodal SGS tensor 5149c678832SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 5150c373b74SJames Wright PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc)); 5169c678832SJames Wright 5179c678832SJames Wright // -- Compute contribution of the SGS stress 518a7dac1d5SJames Wright PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 5199c678832SJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 5209c678832SJames Wright 5219c678832SJames Wright // -- Return local SGS vector 522a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 5239c678832SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 524*8340219bSJames Wright PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 525ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 526d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 5279c678832SJames Wright } 5289c678832SJames Wright 52962b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN 530cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 53162b7942eSJames Wright PetscFunctionBeginUser; 53262b7942eSJames Wright for (PetscInt i = 0; i < N; i++) { 53362b7942eSJames Wright for (PetscInt j = 0; j < M; j++) { 53462b7942eSJames Wright B[j * N + i] = A[i * M + j]; 53562b7942eSJames Wright } 53662b7942eSJames Wright } 537d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 53862b7942eSJames Wright } 53962b7942eSJames Wright 54062b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct 541ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) { 542ad494f68SJames Wright SgsDDContext sgsdd_ctx; 54362b7942eSJames Wright PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 54462b7942eSJames Wright char file_path[PETSC_MAX_PATH_LEN]; 54562b7942eSJames Wright PetscScalar *temp; 54662b7942eSJames Wright 54762b7942eSJames Wright PetscFunctionBeginUser; 54862b7942eSJames Wright { 549ad494f68SJames Wright SgsDDContext sgsdd_temp; 55062b7942eSJames Wright PetscCall(PetscNew(&sgsdd_temp)); 55162b7942eSJames Wright *sgsdd_temp = **psgsdd_ctx; 55262b7942eSJames Wright sgsdd_temp->offsets.bias1 = 0; 55362b7942eSJames Wright sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 55462b7942eSJames Wright sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 55562b7942eSJames Wright sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 55662b7942eSJames Wright sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 55762b7942eSJames Wright PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 55862b7942eSJames Wright sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 55962b7942eSJames Wright PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 56062b7942eSJames Wright *sgsdd_ctx = *sgsdd_temp; 56162b7942eSJames Wright PetscCall(PetscFree(sgsdd_temp)); 56262b7942eSJames Wright } 56362b7942eSJames Wright 56462b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 56542454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 56662b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 56742454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 56862b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 56942454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 57062b7942eSJames Wright 57162b7942eSJames Wright { 57262b7942eSJames Wright PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 57362b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 57442454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 57562b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 57662b7942eSJames Wright PetscCall(PetscFree(temp)); 57762b7942eSJames Wright } 57862b7942eSJames Wright { 57962b7942eSJames Wright PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 58062b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 58142454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 58262b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 58362b7942eSJames Wright PetscCall(PetscFree(temp)); 58462b7942eSJames Wright } 58562b7942eSJames Wright 58662b7942eSJames Wright PetscCall(PetscFree(*psgsdd_ctx)); 58762b7942eSJames Wright *psgsdd_ctx = sgsdd_ctx; 588d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 58962b7942eSJames Wright } 59062b7942eSJames Wright 591e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) { 592ee1455b7SJames Wright PetscReal alpha = 0; 593ad494f68SJames Wright SgsDDContext sgsdd_ctx; 5940c373b74SJames Wright MPI_Comm comm = honee->comm; 595ee1455b7SJames Wright char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 596ad494f68SJames Wright SgsDDSetupData sgs_dd_setup_data; 597ee1455b7SJames Wright NewtonianIdealGasContext gas; 598*8340219bSJames Wright NodalProjectionData grad_velo_proj; 59962b7942eSJames Wright 60006f41313SJames Wright PetscFunctionBeginUser; 601*8340219bSJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &grad_velo_proj)); 602*8340219bSJames Wright PetscCall(PetscObjectContainerCompose((PetscObject)honee, GRAD_VELO_PROJ_KEY, grad_velo_proj, (PetscCtxDestroyFn *)NodalProjectionDataDestroy)); 6039ab09d51SJames Wright 6040c373b74SJames Wright PetscCall(PetscNew(&honee->sgs_dd_data)); 6050c373b74SJames Wright honee->sgs_dd_data->num_comp_inputs = 6; 6060c373b74SJames Wright honee->sgs_dd_data->num_comp_outputs = 6; 60762b7942eSJames Wright 6084c07ec22SJames Wright PetscCall(PetscNew(&sgs_dd_setup_data)); 6094c07ec22SJames Wright 6104b0f6111SJames Wright PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 61162b7942eSJames Wright PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 61262b7942eSJames Wright PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 61362b7942eSJames Wright sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 6144c07ec22SJames Wright PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead")); 6154c07ec22SJames Wright sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED; 6164c07ec22SJames Wright PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations, 6174c07ec22SJames Wright (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation, 6184c07ec22SJames Wright NULL)); 61962b7942eSJames Wright PetscOptionsEnd(); 62062b7942eSJames Wright 621b87d60b3SJames Wright PetscCall(PetscNew(&sgsdd_ctx)); 6224b0f6111SJames Wright sgsdd_ctx->num_layers = 1; 62362b7942eSJames Wright sgsdd_ctx->num_inputs = 6; 62462b7942eSJames Wright sgsdd_ctx->num_outputs = 6; 62562b7942eSJames Wright sgsdd_ctx->num_neurons = 20; 62662b7942eSJames Wright sgsdd_ctx->alpha = alpha; 62762b7942eSJames Wright 628ad494f68SJames Wright PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 62962b7942eSJames Wright 630ee1455b7SJames Wright // -- Create DM for storing SGS tensor at nodes 6310c373b74SJames Wright PetscCall( 6320c373b74SJames Wright SgsDDCreateDM(honee->dm, &honee->sgs_dd_data->dm_sgs, honee->app_ctx->degree, honee->app_ctx->q_extra, &honee->sgs_dd_data->num_comp_sgs)); 633ee1455b7SJames Wright 634e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas)); 635ee1455b7SJames Wright sgsdd_ctx->gas = *gas; 636e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 6370c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 638b4c37c5cSJames Wright PetscCallCeed(ceed, 639b4c37c5cSJames Wright CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 640b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 641ee1455b7SJames Wright 642e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx)); 64340816385SJames Wright 644c38c977aSJames Wright // -- Compute and store anisotropy tensor 645e3663b90SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed)); 646c38c977aSJames Wright 647ee1455b7SJames Wright // -- Create Nodal Evaluation Operator 6484c07ec22SJames Wright switch (sgs_dd_setup_data->sgs_dd_model_implementation) { 6494c07ec22SJames Wright case SGS_MODEL_DD_FUSED: 650e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data)); 6514c07ec22SJames Wright break; 6524c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_CEED: 6534c07ec22SJames Wright case SGS_MODEL_DD_SEQENTIAL_TORCH: 654e3663b90SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data)); 6554c07ec22SJames Wright break; 6564c07ec22SJames Wright } 657ee1455b7SJames Wright 6589c678832SJames Wright // -- Create Operator to evalutate residual of SGS stress 659e3663b90SJames Wright PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data)); 6609c678832SJames Wright 661ad494f68SJames Wright PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data)); 662d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 66362b7942eSJames Wright } 664ee1455b7SJames Wright 66542454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) { 666ee1455b7SJames Wright PetscFunctionBeginUser; 667d949ddfcSJames Wright if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS); 668b4c37c5cSJames Wright Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed; 669ee1455b7SJames Wright 670b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed)); 671b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed)); 672ee1455b7SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx)); 67341edf198SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx)); 674b87d60b3SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx)); 675b87d60b3SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx)); 676ee1455b7SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_sgs)); 677b87d60b3SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs)); 678b87d60b3SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs)); 679b87d60b3SJames 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)); 680ee1455b7SJames Wright PetscCall(PetscFree(sgs_dd_data)); 681d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 682ee1455b7SJames Wright } 683