xref: /honee/problems/sgs_dd_model.c (revision ea615d4cc464aa6ad650c06fae6d120cc2465bc4)
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));
30519781aeSJames Wright   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
31d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32c38c977aSJames Wright }
33c38c977aSJames Wright 
34ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
35ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
36ee1455b7SJames Wright   PetscSection section;
37ee1455b7SJames Wright 
38ee1455b7SJames Wright   PetscFunctionBeginUser;
39ee1455b7SJames Wright   *num_components = 6;
40ee1455b7SJames Wright 
41ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
420dee9b8eSJames Wright   PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE));
43ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
44ee1455b7SJames Wright 
45da4ca0cfSJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
46ee1455b7SJames Wright 
47ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
48ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
49ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
50ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
51ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
52ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
53ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
54ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
55d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
56ee1455b7SJames Wright };
57ee1455b7SJames Wright 
58ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method
590c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
600c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
61cceb3143SJames Wright   PetscMemType q_mem_type;
62cceb3143SJames Wright 
63cceb3143SJames Wright   PetscFunctionBeginUser;
640c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
65cceb3143SJames Wright 
66cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
67cceb3143SJames Wright 
680c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
69cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
70cceb3143SJames Wright }
71cceb3143SJames Wright 
72b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
73e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
740c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
755930f037SJames Wright   CeedQFunction       qf_sgs_dd_nodal;
765930f037SJames Wright   CeedOperator        op_sgs_dd_nodal;
7715c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
78defe8520SJames Wright   PetscInt            dim;
795930f037SJames Wright   CeedVector          inv_multiplicity;
80ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
8115c18037SJames Wright   DMLabel             domain_label = NULL;
8215c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
83ee1455b7SJames Wright 
84ee1455b7SJames Wright   PetscFunctionBeginUser;
850c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
86e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
87e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
88b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
89ee1455b7SJames Wright 
90ee1455b7SJames Wright   {  // Get velocity gradient information
91ee1455b7SJames Wright     CeedOperatorField op_field;
920c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
93b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
94b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
95ee1455b7SJames Wright   }
9615c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
97b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
98ee1455b7SJames Wright 
995930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
1005930f037SJames Wright                                    &inv_multiplicity));
101ee1455b7SJames Wright 
102ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
1030c373b74SJames Wright   switch (honee->phys->state_var) {
104ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
105ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
106ee1455b7SJames Wright       break;
107ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
108ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
109ee1455b7SJames Wright       break;
1109b103f75SJames Wright     case STATEVAR_ENTROPY:
1119b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
1129b103f75SJames Wright       break;
113ee1455b7SJames Wright   }
114ee1455b7SJames Wright 
115ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
116ee1455b7SJames Wright   CeedBasis basis_x_to_q;
117e3663b90SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q));
118ee1455b7SJames Wright 
119b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
120b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
121b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
122b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
123b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
124b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
125b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
126ee1455b7SJames Wright 
127b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
128e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
129e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord));
13058e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
13158e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
132b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
13358e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
13458e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
135ee1455b7SJames Wright 
1360c373b74SJames 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,
1374b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
138ee1455b7SJames Wright 
139ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
140ad494f68SJames Wright   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
141ee1455b7SJames Wright 
142b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
143b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
144b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
145fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
146b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
147b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
148d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
149ee1455b7SJames Wright }
150ee1455b7SJames Wright 
1514c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation
1524c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1534c07ec22SJames Wright                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
154b87d60b3SJames Wright                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
155b87d60b3SJames Wright                                                                 void **ctx) {
156b87d60b3SJames Wright   CeedQFunction         qf_sgs_dd_inference;
157b87d60b3SJames Wright   CeedOperator          op_sgs_dd_inference;
158b87d60b3SJames Wright   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
159b87d60b3SJames Wright 
160b87d60b3SJames Wright   PetscFunctionBeginUser;
161b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
162b87d60b3SJames Wright                                                   &qf_sgs_dd_inference));
163b87d60b3SJames Wright 
164b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
165b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
166b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
167b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
168b87d60b3SJames Wright 
169b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
170b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
171b87d60b3SJames Wright   PetscCallCeed(ceed,
172b87d60b3SJames Wright                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
173b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
174b87d60b3SJames Wright 
175b87d60b3SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
176b87d60b3SJames Wright                                        op_context));
177b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
178b87d60b3SJames Wright 
179b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
180b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
181b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
182b87d60b3SJames Wright }
183b87d60b3SJames Wright 
1844c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation
1854c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
186b87d60b3SJames Wright   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
187b87d60b3SJames Wright 
188b87d60b3SJames Wright   PetscFunctionBeginUser;
189*ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
190*ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
191b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeBegin());
192b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
193b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeEnd());
194*ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
195*ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
196b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
197b87d60b3SJames Wright }
198b87d60b3SJames Wright 
1994c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch
2004c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
2014c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
2024c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
2034c07ec22SJames Wright                                                                  void **ctx) {
2044c07ec22SJames Wright   const char     *ceed_resource;
2054c07ec22SJames Wright   char            model_path[PETSC_MAX_PATH_LEN] = "";
2064c07ec22SJames Wright   TorchDeviceType model_device_type;
2074c07ec22SJames Wright 
2084c07ec22SJames Wright   PetscFunctionBeginUser;
2094c07ec22SJames Wright   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
2104c07ec22SJames Wright   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
2114c07ec22SJames Wright   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
2127ffa0ff8SJames Wright   // On-device XPU is not working reliably currently, default to CPU inference evaluation
2137ffa0ff8SJames Wright   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
2144c07ec22SJames Wright   else model_device_type = TORCH_DEVICE_CPU;
2154c07ec22SJames Wright   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
2164c07ec22SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
2174c07ec22SJames Wright 
2184c07ec22SJames Wright   PetscCall(LoadModel_Torch(model_path, model_device_type));
2194c07ec22SJames Wright 
2204c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2214c07ec22SJames Wright }
2224c07ec22SJames Wright 
2234c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch
2244c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
2254c07ec22SJames Wright   static PetscBool run_through = PETSC_FALSE;
2264c07ec22SJames Wright   PetscFunctionBeginUser;
2274c07ec22SJames Wright   if (!run_through) {
2284c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
2294c07ec22SJames Wright   }
2304c07ec22SJames Wright   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
2314c07ec22SJames Wright   if (!run_through) {
2324c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
2334c07ec22SJames Wright     run_through = PETSC_TRUE;
2344c07ec22SJames Wright   }
2354c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2364c07ec22SJames Wright }
2374c07ec22SJames Wright 
238b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method
2390c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
2400c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
241b87d60b3SJames Wright   PetscMemType q_mem_type;
242b87d60b3SJames Wright   Vec          DD_Inputs_loc, DD_Outputs_loc;
243b87d60b3SJames Wright 
244b87d60b3SJames Wright   PetscFunctionBeginUser;
245b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
246b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
2470c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
248b87d60b3SJames Wright 
249b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
250b87d60b3SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
251b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
252b87d60b3SJames Wright 
2530c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
254b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
255b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
256b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
257b87d60b3SJames Wright }
258b87d60b3SJames Wright 
259b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
260e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
2610c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
262b87d60b3SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
263b87d60b3SJames Wright   PetscInt            dim;
264b87d60b3SJames Wright   CeedVector          inv_multiplicity, eigvec;
265b87d60b3SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
266b87d60b3SJames Wright       elem_restr_dd_outputs;
267b87d60b3SJames Wright   DMLabel  domain_label = NULL;
268b87d60b3SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
269b87d60b3SJames Wright 
270b87d60b3SJames Wright   PetscFunctionBeginUser;
271b87d60b3SJames Wright   {  // Create DMs for data-driven input and output values
272b87d60b3SJames Wright     PetscSection section;
273b87d60b3SJames Wright     PetscInt     degree, q_extra;
274b87d60b3SJames Wright     {  // Get degree and number of quadrature points from dm_sgs
275b87d60b3SJames Wright       PetscFE         fe;
276b87d60b3SJames Wright       PetscSpace      basis;
277b87d60b3SJames Wright       PetscQuadrature quadrature;
278b87d60b3SJames Wright       PetscInt        num_qpnts;
279b87d60b3SJames Wright       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
280b87d60b3SJames Wright       PetscCall(PetscFEGetBasisSpace(fe, &basis));
281b87d60b3SJames Wright       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
282b87d60b3SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
283b87d60b3SJames Wright       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
284b87d60b3SJames Wright       q_extra = degree - num_qpnts;
285b87d60b3SJames Wright     }
286b87d60b3SJames Wright 
287b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
2880dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE));
289b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
290b87d60b3SJames 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));
291b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
292b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
293b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
294b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
295b87d60b3SJames Wright 
296b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
297b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
298b87d60b3SJames Wright     }
299b87d60b3SJames Wright 
300b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
3010dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE));
302b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
303b87d60b3SJames 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));
304b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
305b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
306b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
307b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
308b87d60b3SJames Wright 
309b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
310b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
311b87d60b3SJames Wright     }
312b87d60b3SJames Wright   }
313b87d60b3SJames Wright 
3140c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
315e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
316e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
317b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
318b87d60b3SJames Wright 
319b87d60b3SJames Wright   {  // Get velocity gradient information
320b87d60b3SJames Wright     CeedOperatorField op_field;
3210c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
322b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
323b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
324b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
325b87d60b3SJames Wright   }
326b87d60b3SJames Wright 
327b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
328b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
329b87d60b3SJames Wright   PetscCall(
330b87d60b3SJames Wright       DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec));
331b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
332b87d60b3SJames Wright 
333b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs));
334b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs));
335b87d60b3SJames Wright 
3365930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
3375930f037SJames Wright                                    &inv_multiplicity));
338b87d60b3SJames Wright 
339b87d60b3SJames Wright   {  // Create operator for data-driven input evaluation
340b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_inputs;
341b87d60b3SJames Wright     CeedOperator  op_sgs_dd_inputs;
342b87d60b3SJames Wright 
3430c373b74SJames Wright     switch (honee->phys->state_var) {
344b87d60b3SJames Wright       case STATEVAR_PRIMITIVE:
345b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
346b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
347b87d60b3SJames Wright         break;
348b87d60b3SJames Wright       case STATEVAR_CONSERVATIVE:
349b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
350b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
351b87d60b3SJames Wright         break;
3529b103f75SJames Wright       case STATEVAR_ENTROPY:
3539b103f75SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
3549b103f75SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
3559b103f75SJames Wright         break;
356b87d60b3SJames Wright     }
357b87d60b3SJames Wright 
358b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
359b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
360b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
361b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
362b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
363b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
364b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
365b87d60b3SJames Wright 
366b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
367e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
368b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
369b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
370b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
371b87d60b3SJames Wright     PetscCallCeed(ceed,
372b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
373b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
374b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
375b87d60b3SJames Wright 
3760c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(honee->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
377b87d60b3SJames Wright                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
378b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
379b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
380b87d60b3SJames Wright   }
381b87d60b3SJames Wright 
382b87d60b3SJames Wright   {  // Create operator for data-driven output handling
383b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_outputs;
384b87d60b3SJames Wright     CeedOperator  op_sgs_dd_outputs;
385b87d60b3SJames Wright 
386b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
387b87d60b3SJames Wright                                                     &qf_sgs_dd_outputs));
388b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
389b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
390b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
391b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
392b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
393b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
394b87d60b3SJames Wright 
395b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
396b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
397b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
398b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
399b87d60b3SJames Wright     PetscCallCeed(ceed,
400b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
401b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
402b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
403b87d60b3SJames Wright 
404b87d60b3SJames 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,
405b87d60b3SJames Wright                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
406b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
407b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
408b87d60b3SJames Wright   }
409b87d60b3SJames Wright 
410b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
4114c07ec22SJames Wright 
4124c07ec22SJames Wright   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
4134c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
4144c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
415b87d60b3SJames Wright                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4164c07ec22SJames Wright   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
4174c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
4184c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
4194c07ec22SJames Wright                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4204c07ec22SJames Wright   }
421b87d60b3SJames Wright 
422b87d60b3SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
423b87d60b3SJames Wright 
424b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
425b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
426b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
427b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
428b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
429b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
430fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
431b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
432b87d60b3SJames Wright }
433b87d60b3SJames Wright 
4349c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
435e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
4360c373b74SJames Wright   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
437be29160dSJames Wright   CeedInt             num_comp_q, q_data_size, num_comp_x;
438defe8520SJames Wright   PetscInt            dim;
4399c678832SJames Wright   CeedQFunction       qf_sgs_apply;
4409c678832SJames Wright   CeedOperator        op_sgs_apply;
4419c678832SJames Wright   CeedBasis           basis_sgs;
442be29160dSJames Wright   CeedVector          q_data;
443be29160dSJames Wright   CeedElemRestriction elem_restr_qd;
4449c678832SJames Wright 
4459c678832SJames Wright   PetscFunctionBeginUser;
4460c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
447e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
448e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
4499c678832SJames Wright 
450be29160dSJames Wright   {
451be29160dSJames Wright     DMLabel  domain_label = NULL;
452be29160dSJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
453be29160dSJames Wright 
454be29160dSJames Wright     PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs));
455e3663b90SJames 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,
456e3663b90SJames Wright                        &q_data, &q_data_size));
457be29160dSJames Wright   }
4589c678832SJames Wright 
4590c373b74SJames Wright   switch (honee->phys->state_var) {
4609c678832SJames Wright     case STATEVAR_PRIMITIVE:
46142454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
4629c678832SJames Wright       break;
4639c678832SJames Wright     case STATEVAR_CONSERVATIVE:
46442454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
4659c678832SJames Wright       break;
4669b103f75SJames Wright     case STATEVAR_ENTROPY:
4679b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
4689b103f75SJames Wright       break;
4699c678832SJames Wright   }
4709c678832SJames Wright 
47140816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
472b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
473be29160dSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE));
474b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
475b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
4769c678832SJames Wright 
477b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
478e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
479be29160dSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
480b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
481e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
4829c678832SJames Wright 
4839c678832SJames Wright   PetscCall(
4840c373b74SJames 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));
4859c678832SJames Wright 
486be29160dSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
487be29160dSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
48887edb941SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs));
489b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
490b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
491d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4929c678832SJames Wright }
4939c678832SJames Wright 
4949c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
4950c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) {
4960c373b74SJames Wright   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
4979c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
498cceb3143SJames Wright   PetscMemType sgs_nodal_mem_type;
4999c678832SJames Wright 
5009c678832SJames Wright   PetscFunctionBeginUser;
501*ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL));
5020c373b74SJames Wright   PetscCall(DMGetGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
5030c373b74SJames Wright   PetscCall(VelocityGradientProjectionApply(honee->grad_velo_proj, Q_loc, VelocityGradient));
5049c678832SJames Wright 
5059c678832SJames Wright   // -- Compute Nodal SGS tensor
5069c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5070c373b74SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc));
5089c678832SJames Wright 
5099c678832SJames Wright   // -- Compute contribution of the SGS stress
510a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
5119c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
5129c678832SJames Wright 
5139c678832SJames Wright   // -- Return local SGS vector
514a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
5159c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5160c373b74SJames Wright   PetscCall(DMRestoreGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
517*ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL));
518d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5199c678832SJames Wright }
5209c678832SJames Wright 
52162b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
522cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
52362b7942eSJames Wright   PetscFunctionBeginUser;
52462b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
52562b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
52662b7942eSJames Wright       B[j * N + i] = A[i * M + j];
52762b7942eSJames Wright     }
52862b7942eSJames Wright   }
529d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53062b7942eSJames Wright }
53162b7942eSJames Wright 
53262b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
533ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
534ad494f68SJames Wright   SgsDDContext sgsdd_ctx;
53562b7942eSJames Wright   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
53662b7942eSJames Wright   char         file_path[PETSC_MAX_PATH_LEN];
53762b7942eSJames Wright   PetscScalar *temp;
53862b7942eSJames Wright 
53962b7942eSJames Wright   PetscFunctionBeginUser;
54062b7942eSJames Wright   {
541ad494f68SJames Wright     SgsDDContext sgsdd_temp;
54262b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
54362b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
54462b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
54562b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
54662b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
54762b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
54862b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
54962b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
55062b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
55162b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
55262b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
55362b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
55462b7942eSJames Wright   }
55562b7942eSJames Wright 
55662b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
55742454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
55862b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
55942454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
56062b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
56142454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
56262b7942eSJames Wright 
56362b7942eSJames Wright   {
56462b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
56562b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
56642454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
56762b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
56862b7942eSJames Wright     PetscCall(PetscFree(temp));
56962b7942eSJames Wright   }
57062b7942eSJames Wright   {
57162b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
57262b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
57342454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
57462b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
57562b7942eSJames Wright     PetscCall(PetscFree(temp));
57662b7942eSJames Wright   }
57762b7942eSJames Wright 
57862b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
57962b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
580d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
58162b7942eSJames Wright }
58262b7942eSJames Wright 
583e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) {
584ee1455b7SJames Wright   PetscReal                alpha = 0;
585ad494f68SJames Wright   SgsDDContext             sgsdd_ctx;
5860c373b74SJames Wright   MPI_Comm                 comm                           = honee->comm;
587ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
588ad494f68SJames Wright   SgsDDSetupData           sgs_dd_setup_data;
589ee1455b7SJames Wright   NewtonianIdealGasContext gas;
59062b7942eSJames Wright 
59106f41313SJames Wright   PetscFunctionBeginUser;
592e3663b90SJames Wright   PetscCall(
593e3663b90SJames Wright       VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &honee->grad_velo_proj));
5949ab09d51SJames Wright 
5950c373b74SJames Wright   PetscCall(PetscNew(&honee->sgs_dd_data));
5960c373b74SJames Wright   honee->sgs_dd_data->num_comp_inputs  = 6;
5970c373b74SJames Wright   honee->sgs_dd_data->num_comp_outputs = 6;
59862b7942eSJames Wright 
5994c07ec22SJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
6004c07ec22SJames Wright 
6014b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
60262b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
60362b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
60462b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
6054c07ec22SJames Wright   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
6064c07ec22SJames Wright   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
6074c07ec22SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
6084c07ec22SJames Wright                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
6094c07ec22SJames Wright                              NULL));
61062b7942eSJames Wright   PetscOptionsEnd();
61162b7942eSJames Wright 
612b87d60b3SJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
6134b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
61462b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
61562b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
61662b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
61762b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
61862b7942eSJames Wright 
619ad494f68SJames Wright   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
62062b7942eSJames Wright 
621ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
6220c373b74SJames Wright   PetscCall(
6230c373b74SJames 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));
624ee1455b7SJames Wright 
625e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
626ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
627e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
6280c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
629b4c37c5cSJames Wright   PetscCallCeed(ceed,
630b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
631b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
632ee1455b7SJames Wright 
633e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx));
63440816385SJames Wright 
635c38c977aSJames Wright   // -- Compute and store anisotropy tensor
636e3663b90SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed));
637c38c977aSJames Wright 
638ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
6394c07ec22SJames Wright   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
6404c07ec22SJames Wright     case SGS_MODEL_DD_FUSED:
641e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data));
6424c07ec22SJames Wright       break;
6434c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_CEED:
6444c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_TORCH:
645e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data));
6464c07ec22SJames Wright       break;
6474c07ec22SJames Wright   }
648ee1455b7SJames Wright 
6499c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
650e3663b90SJames Wright   PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data));
6519c678832SJames Wright 
652ad494f68SJames Wright   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
653d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65462b7942eSJames Wright }
655ee1455b7SJames Wright 
65642454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
657ee1455b7SJames Wright   PetscFunctionBeginUser;
658d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
659b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
660ee1455b7SJames Wright 
661b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
662b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed));
663ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
66441edf198SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
665b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx));
666b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx));
667ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
668b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs));
669b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs));
670b87d60b3SJames 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));
671ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
672d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
673ee1455b7SJames Wright }
674