xref: /honee/problems/sgs_dd_model.c (revision fff85bd3a0c82b434cd3f01285a22a362dfa8070)
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 
18ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
19b4c37c5cSJames Wright   Ceed ceed;
20ad494f68SJames Wright 
21c38c977aSJames Wright   PetscFunctionBeginUser;
22b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
23ad494f68SJames Wright 
24b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
25b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
26b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
27b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
2840816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
29c38c977aSJames Wright   PetscCall(PetscFree(sgs_dd_setup_data));
30d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
31c38c977aSJames Wright }
32c38c977aSJames Wright 
33ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
34ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
35ee1455b7SJames Wright   PetscSection section;
36ee1455b7SJames Wright 
37ee1455b7SJames Wright   PetscFunctionBeginUser;
38ee1455b7SJames Wright   *num_components = 6;
39ee1455b7SJames Wright 
40ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
410dee9b8eSJames Wright   PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE));
42ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
43ee1455b7SJames Wright 
44da4ca0cfSJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
45ee1455b7SJames Wright 
46ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
47ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
48ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
49ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
50ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
51ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
52ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
53ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
54d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
55ee1455b7SJames Wright };
56ee1455b7SJames Wright 
57ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method
580c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
590c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
60cceb3143SJames Wright   PetscMemType q_mem_type;
61cceb3143SJames Wright 
62cceb3143SJames Wright   PetscFunctionBeginUser;
630c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
64cceb3143SJames Wright 
65cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
66cceb3143SJames Wright 
670c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
68cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
69cceb3143SJames Wright }
70cceb3143SJames Wright 
71b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
72e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
730c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
745930f037SJames Wright   CeedQFunction       qf_sgs_dd_nodal;
755930f037SJames Wright   CeedOperator        op_sgs_dd_nodal;
7615c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
77defe8520SJames Wright   PetscInt            dim;
785930f037SJames Wright   CeedVector          inv_multiplicity;
79ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
8015c18037SJames Wright   DMLabel             domain_label = NULL;
8115c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
82ee1455b7SJames Wright 
83ee1455b7SJames Wright   PetscFunctionBeginUser;
840c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
85e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
86e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
87b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
88ee1455b7SJames Wright 
89ee1455b7SJames Wright   {  // Get velocity gradient information
90ee1455b7SJames Wright     CeedOperatorField op_field;
910c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
92b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
93b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
94ee1455b7SJames Wright   }
9515c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
96b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
97ee1455b7SJames Wright 
985930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
995930f037SJames Wright                                    &inv_multiplicity));
100ee1455b7SJames Wright 
101ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
1020c373b74SJames Wright   switch (honee->phys->state_var) {
103ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
104ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
105ee1455b7SJames Wright       break;
106ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
107ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
108ee1455b7SJames Wright       break;
1099b103f75SJames Wright     case STATEVAR_ENTROPY:
1109b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
1119b103f75SJames Wright       break;
112ee1455b7SJames Wright   }
113ee1455b7SJames Wright 
114ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
115ee1455b7SJames Wright   CeedBasis basis_x_to_q;
116e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q));
117ee1455b7SJames Wright 
118b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
119b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
120b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
121b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
122b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
123b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
124b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
125ee1455b7SJames Wright 
126b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
127e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
128e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord));
12958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
13058e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
131b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
13258e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
13358e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
134ee1455b7SJames Wright 
1350c373b74SJames Wright   PetscCall(OperatorApplyContextCreate(honee->grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, NULL, sgs_dd_data->sgs_nodal_ceed, NULL,
1364b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
137ee1455b7SJames Wright 
138ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
139ad494f68SJames Wright   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
140ee1455b7SJames Wright 
141b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
142b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
143b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
144*fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
145b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
146b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
147d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
148ee1455b7SJames Wright }
149ee1455b7SJames Wright 
1504c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation
1514c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1524c07ec22SJames Wright                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
153b87d60b3SJames Wright                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
154b87d60b3SJames Wright                                                                 void **ctx) {
155b87d60b3SJames Wright   CeedQFunction         qf_sgs_dd_inference;
156b87d60b3SJames Wright   CeedOperator          op_sgs_dd_inference;
157b87d60b3SJames Wright   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
158b87d60b3SJames Wright 
159b87d60b3SJames Wright   PetscFunctionBeginUser;
160b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
161b87d60b3SJames Wright                                                   &qf_sgs_dd_inference));
162b87d60b3SJames Wright 
163b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
164b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
165b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
166b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
167b87d60b3SJames Wright 
168b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
169b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
170b87d60b3SJames Wright   PetscCallCeed(ceed,
171b87d60b3SJames Wright                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
172b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
173b87d60b3SJames Wright 
174b87d60b3SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
175b87d60b3SJames Wright                                        op_context));
176b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
177b87d60b3SJames Wright 
178b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
179b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
180b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
181b87d60b3SJames Wright }
182b87d60b3SJames Wright 
1834c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation
1844c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
185b87d60b3SJames Wright   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
186b87d60b3SJames Wright 
187b87d60b3SJames Wright   PetscFunctionBeginUser;
188b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
189b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
190b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeBegin());
191b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
192b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeEnd());
193b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
194b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
195b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
196b87d60b3SJames Wright }
197b87d60b3SJames Wright 
1984c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch
1994c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
2004c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
2014c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
2024c07ec22SJames Wright                                                                  void **ctx) {
2034c07ec22SJames Wright   const char     *ceed_resource;
2044c07ec22SJames Wright   char            model_path[PETSC_MAX_PATH_LEN] = "";
2054c07ec22SJames Wright   TorchDeviceType model_device_type;
2064c07ec22SJames Wright 
2074c07ec22SJames Wright   PetscFunctionBeginUser;
2084c07ec22SJames Wright   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
2094c07ec22SJames Wright   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
2104c07ec22SJames Wright   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
2117ffa0ff8SJames Wright   // On-device XPU is not working reliably currently, default to CPU inference evaluation
2127ffa0ff8SJames Wright   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
2134c07ec22SJames Wright   else model_device_type = TORCH_DEVICE_CPU;
2144c07ec22SJames Wright   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
2154c07ec22SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
2164c07ec22SJames Wright 
2174c07ec22SJames Wright   PetscCall(LoadModel_Torch(model_path, model_device_type));
2184c07ec22SJames Wright 
2194c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2204c07ec22SJames Wright }
2214c07ec22SJames Wright 
2224c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch
2234c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
2244c07ec22SJames Wright   static PetscBool run_through = PETSC_FALSE;
2254c07ec22SJames Wright   PetscFunctionBeginUser;
2264c07ec22SJames Wright   if (!run_through) {
2274c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
2284c07ec22SJames Wright   }
2294c07ec22SJames Wright   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
2304c07ec22SJames Wright   if (!run_through) {
2314c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
2324c07ec22SJames Wright     run_through = PETSC_TRUE;
2334c07ec22SJames Wright   }
2344c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2354c07ec22SJames Wright }
2364c07ec22SJames Wright 
237b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method
2380c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
2390c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
240b87d60b3SJames Wright   PetscMemType q_mem_type;
241b87d60b3SJames Wright   Vec          DD_Inputs_loc, DD_Outputs_loc;
242b87d60b3SJames Wright 
243b87d60b3SJames Wright   PetscFunctionBeginUser;
244b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
245b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
2460c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
247b87d60b3SJames Wright 
248b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
249b87d60b3SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
250b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
251b87d60b3SJames Wright 
2520c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
253b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
254b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
255b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
256b87d60b3SJames Wright }
257b87d60b3SJames Wright 
258b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
259e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
2600c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
261b87d60b3SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
262b87d60b3SJames Wright   PetscInt            dim;
263b87d60b3SJames Wright   CeedVector          inv_multiplicity, eigvec;
264b87d60b3SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
265b87d60b3SJames Wright       elem_restr_dd_outputs;
266b87d60b3SJames Wright   DMLabel  domain_label = NULL;
267b87d60b3SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
268b87d60b3SJames Wright 
269b87d60b3SJames Wright   PetscFunctionBeginUser;
270b87d60b3SJames Wright   {  // Create DMs for data-driven input and output values
271b87d60b3SJames Wright     PetscSection section;
272b87d60b3SJames Wright     PetscInt     degree, q_extra;
273b87d60b3SJames Wright     {  // Get degree and number of quadrature points from dm_sgs
274b87d60b3SJames Wright       PetscFE         fe;
275b87d60b3SJames Wright       PetscSpace      basis;
276b87d60b3SJames Wright       PetscQuadrature quadrature;
277b87d60b3SJames Wright       PetscInt        num_qpnts;
278b87d60b3SJames Wright       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
279b87d60b3SJames Wright       PetscCall(PetscFEGetBasisSpace(fe, &basis));
280b87d60b3SJames Wright       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
281b87d60b3SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
282b87d60b3SJames Wright       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
283b87d60b3SJames Wright       q_extra = degree - num_qpnts;
284b87d60b3SJames Wright     }
285b87d60b3SJames Wright 
286b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
2870dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE));
288b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
289b87d60b3SJames 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));
290b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
291b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
292b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
293b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
294b87d60b3SJames Wright 
295b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
296b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
297b87d60b3SJames Wright     }
298b87d60b3SJames Wright 
299b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
3000dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE));
301b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
302b87d60b3SJames 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));
303b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
304b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
305b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
306b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
307b87d60b3SJames Wright 
308b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
309b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
310b87d60b3SJames Wright     }
311b87d60b3SJames Wright   }
312b87d60b3SJames Wright 
3130c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
314e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
315e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
316b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
317b87d60b3SJames Wright 
318b87d60b3SJames Wright   {  // Get velocity gradient information
319b87d60b3SJames Wright     CeedOperatorField op_field;
3200c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
321b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
322b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
323b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
324b87d60b3SJames Wright   }
325b87d60b3SJames Wright 
326b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
327b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
328b87d60b3SJames Wright   PetscCall(
329b87d60b3SJames Wright       DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec));
330b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
331b87d60b3SJames Wright 
332b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs));
333b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs));
334b87d60b3SJames Wright 
3355930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
3365930f037SJames Wright                                    &inv_multiplicity));
337b87d60b3SJames Wright 
338b87d60b3SJames Wright   {  // Create operator for data-driven input evaluation
339b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_inputs;
340b87d60b3SJames Wright     CeedOperator  op_sgs_dd_inputs;
341b87d60b3SJames Wright 
3420c373b74SJames Wright     switch (honee->phys->state_var) {
343b87d60b3SJames Wright       case STATEVAR_PRIMITIVE:
344b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
345b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
346b87d60b3SJames Wright         break;
347b87d60b3SJames Wright       case STATEVAR_CONSERVATIVE:
348b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
349b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
350b87d60b3SJames Wright         break;
3519b103f75SJames Wright       case STATEVAR_ENTROPY:
3529b103f75SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
3539b103f75SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
3549b103f75SJames Wright         break;
355b87d60b3SJames Wright     }
356b87d60b3SJames Wright 
357b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
358b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
359b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
360b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
361b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
362b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
363b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
364b87d60b3SJames Wright 
365b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
366e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
367b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
368b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
369b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
370b87d60b3SJames Wright     PetscCallCeed(ceed,
371b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
372b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
373b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
374b87d60b3SJames Wright 
3750c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
376b87d60b3SJames Wright                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
377b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
378b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
379b87d60b3SJames Wright   }
380b87d60b3SJames Wright 
381b87d60b3SJames Wright   {  // Create operator for data-driven output handling
382b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_outputs;
383b87d60b3SJames Wright     CeedOperator  op_sgs_dd_outputs;
384b87d60b3SJames Wright 
385b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
386b87d60b3SJames Wright                                                     &qf_sgs_dd_outputs));
387b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
388b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
389b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
390b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
391b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
392b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
393b87d60b3SJames Wright 
394b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
395b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
396b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
397b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
398b87d60b3SJames Wright     PetscCallCeed(ceed,
399b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
400b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
401b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
402b87d60b3SJames Wright 
403b87d60b3SJames 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,
404b87d60b3SJames Wright                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
405b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
406b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
407b87d60b3SJames Wright   }
408b87d60b3SJames Wright 
409b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
4104c07ec22SJames Wright 
4114c07ec22SJames Wright   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
4124c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
4134c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
414b87d60b3SJames Wright                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4154c07ec22SJames Wright   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
4164c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
4174c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
4184c07ec22SJames Wright                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4194c07ec22SJames Wright   }
420b87d60b3SJames Wright 
421b87d60b3SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
422b87d60b3SJames Wright 
423b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
424b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
425b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
426b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
427b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
428b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
429*fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
430b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
431b87d60b3SJames Wright }
432b87d60b3SJames Wright 
4339c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
434e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
4350c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
436be29160dSJames Wright   CeedInt             num_comp_q, q_data_size, num_comp_x;
437defe8520SJames Wright   PetscInt            dim;
4389c678832SJames Wright   CeedQFunction       qf_sgs_apply;
4399c678832SJames Wright   CeedOperator        op_sgs_apply;
4409c678832SJames Wright   CeedBasis           basis_sgs;
441be29160dSJames Wright   CeedVector          q_data;
442be29160dSJames Wright   CeedElemRestriction elem_restr_qd;
4439c678832SJames Wright 
4449c678832SJames Wright   PetscFunctionBeginUser;
4450c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
446e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
447e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
4489c678832SJames Wright 
449be29160dSJames Wright   {
450be29160dSJames Wright     DMLabel  domain_label = NULL;
451be29160dSJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
452be29160dSJames Wright 
453be29160dSJames Wright     PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs));
454e3663b90SJames 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,
455e3663b90SJames Wright                        &q_data, &q_data_size));
456be29160dSJames Wright   }
4579c678832SJames Wright 
4580c373b74SJames Wright   switch (honee->phys->state_var) {
4599c678832SJames Wright     case STATEVAR_PRIMITIVE:
46042454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
4619c678832SJames Wright       break;
4629c678832SJames Wright     case STATEVAR_CONSERVATIVE:
46342454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
4649c678832SJames Wright       break;
4659b103f75SJames Wright     case STATEVAR_ENTROPY:
4669b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
4679b103f75SJames Wright       break;
4689c678832SJames Wright   }
4699c678832SJames Wright 
47040816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
471b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
472be29160dSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE));
473b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
474b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
4759c678832SJames Wright 
476b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
477e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
478be29160dSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
479b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
480e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
4819c678832SJames Wright 
4829c678832SJames Wright   PetscCall(
4830c373b74SJames 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));
4849c678832SJames Wright 
485be29160dSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
486be29160dSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
48787edb941SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs));
488b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
489b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
490d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4919c678832SJames Wright }
4929c678832SJames Wright 
4939c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
4940c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) {
4950c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
4969c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
497cceb3143SJames Wright   PetscMemType sgs_nodal_mem_type;
4989c678832SJames Wright 
4999c678832SJames Wright   PetscFunctionBeginUser;
500b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
5010c373b74SJames Wright   PetscCall(DMGetGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
5020c373b74SJames Wright   PetscCall(VelocityGradientProjectionApply(honee->grad_velo_proj, Q_loc, VelocityGradient));
5039c678832SJames Wright 
5049c678832SJames Wright   // -- Compute Nodal SGS tensor
5059c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5060c373b74SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc));
5079c678832SJames Wright 
5089c678832SJames Wright   // -- Compute contribution of the SGS stress
509a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
5109c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
5119c678832SJames Wright 
5129c678832SJames Wright   // -- Return local SGS vector
513a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
5149c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5150c373b74SJames Wright   PetscCall(DMRestoreGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
516b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
517d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5189c678832SJames Wright }
5199c678832SJames Wright 
52062b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
521cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
52262b7942eSJames Wright   PetscFunctionBeginUser;
52362b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
52462b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
52562b7942eSJames Wright       B[j * N + i] = A[i * M + j];
52662b7942eSJames Wright     }
52762b7942eSJames Wright   }
528d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
52962b7942eSJames Wright }
53062b7942eSJames Wright 
53162b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
532ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
533ad494f68SJames Wright   SgsDDContext sgsdd_ctx;
53462b7942eSJames Wright   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
53562b7942eSJames Wright   char         file_path[PETSC_MAX_PATH_LEN];
53662b7942eSJames Wright   PetscScalar *temp;
53762b7942eSJames Wright 
53862b7942eSJames Wright   PetscFunctionBeginUser;
53962b7942eSJames Wright   {
540ad494f68SJames Wright     SgsDDContext sgsdd_temp;
54162b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
54262b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
54362b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
54462b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
54562b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
54662b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
54762b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
54862b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
54962b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
55062b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
55162b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
55262b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
55362b7942eSJames Wright   }
55462b7942eSJames Wright 
55562b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
55642454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
55762b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
55842454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
55962b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
56042454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
56162b7942eSJames Wright 
56262b7942eSJames Wright   {
56362b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
56462b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
56542454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
56662b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
56762b7942eSJames Wright     PetscCall(PetscFree(temp));
56862b7942eSJames Wright   }
56962b7942eSJames Wright   {
57062b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
57162b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
57242454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
57362b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
57462b7942eSJames Wright     PetscCall(PetscFree(temp));
57562b7942eSJames Wright   }
57662b7942eSJames Wright 
57762b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
57862b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
579d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
58062b7942eSJames Wright }
58162b7942eSJames Wright 
582e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) {
583ee1455b7SJames Wright   PetscReal                alpha = 0;
584ad494f68SJames Wright   SgsDDContext             sgsdd_ctx;
5850c373b74SJames Wright   MPI_Comm                 comm                           = honee->comm;
586ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
587ad494f68SJames Wright   SgsDDSetupData           sgs_dd_setup_data;
588ee1455b7SJames Wright   NewtonianIdealGasContext gas;
58962b7942eSJames Wright 
59006f41313SJames Wright   PetscFunctionBeginUser;
591e3663b90SJames Wright   PetscCall(
592e3663b90SJames Wright       VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &honee->grad_velo_proj));
5939ab09d51SJames Wright 
5940c373b74SJames Wright   PetscCall(PetscNew(&honee->sgs_dd_data));
5950c373b74SJames Wright   honee->sgs_dd_data->num_comp_inputs  = 6;
5960c373b74SJames Wright   honee->sgs_dd_data->num_comp_outputs = 6;
59762b7942eSJames Wright 
5984c07ec22SJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
5994c07ec22SJames Wright 
6004b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
60162b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
60262b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
60362b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
6044c07ec22SJames Wright   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
6054c07ec22SJames Wright   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
6064c07ec22SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
6074c07ec22SJames Wright                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
6084c07ec22SJames Wright                              NULL));
60962b7942eSJames Wright   PetscOptionsEnd();
61062b7942eSJames Wright 
611b87d60b3SJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
6124b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
61362b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
61462b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
61562b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
61662b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
61762b7942eSJames Wright 
618ad494f68SJames Wright   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
61962b7942eSJames Wright 
620ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
6210c373b74SJames Wright   PetscCall(
6220c373b74SJames 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));
623ee1455b7SJames Wright 
624e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
625ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
626e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
6270c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
628b4c37c5cSJames Wright   PetscCallCeed(ceed,
629b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
630b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
631ee1455b7SJames Wright 
632e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx));
63340816385SJames Wright 
634c38c977aSJames Wright   // -- Compute and store anisotropy tensor
635e3663b90SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed));
636c38c977aSJames Wright 
637ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
6384c07ec22SJames Wright   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
6394c07ec22SJames Wright     case SGS_MODEL_DD_FUSED:
640e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data));
6414c07ec22SJames Wright       break;
6424c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_CEED:
6434c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_TORCH:
644e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data));
6454c07ec22SJames Wright       break;
6464c07ec22SJames Wright   }
647ee1455b7SJames Wright 
6489c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
649e3663b90SJames Wright   PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data));
6509c678832SJames Wright 
651ad494f68SJames Wright   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
652d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65362b7942eSJames Wright }
654ee1455b7SJames Wright 
65542454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
656ee1455b7SJames Wright   PetscFunctionBeginUser;
657d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
658b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
659ee1455b7SJames Wright 
660b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
661b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed));
662ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
66341edf198SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
664b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx));
665b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx));
666ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
667b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs));
668b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs));
669b87d60b3SJames 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));
670ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
671d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
672ee1455b7SJames Wright }
673