xref: /honee/problems/sgs_dd_model.c (revision 67aa9f91519d5d2e65314c9a1a411fa9d60d3e06)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 #include "../qfunctions/sgs_dd_model.h"
5 
6 #include <petscdmplex.h>
7 
8 #include <navierstokes.h>
9 #include <sgs_model_torch.h>
10 
11 typedef struct {
12   CeedElemRestriction      elem_restr_grid_aniso, elem_restr_sgs;
13   CeedVector               grid_aniso_ceed;
14   CeedQFunctionContext     sgsdd_qfctx, ifunction_qfctx;
15   SGSModelDDImplementation sgs_dd_model_implementation;
16 } *SgsDDSetupData;
17 
18 PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
19   Ceed ceed;
20 
21   PetscFunctionBeginUser;
22   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
23 
24   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
25   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
26   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
27   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
28   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
29   PetscCall(PetscFree(sgs_dd_setup_data));
30   PetscFunctionReturn(PETSC_SUCCESS);
31 }
32 
33 // @brief Create DM for storing subgrid stress at nodes
34 static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
35   PetscSection section;
36 
37   PetscFunctionBeginUser;
38   *num_components = 6;
39 
40   PetscCall(DMClone(dm_source, dm_sgs));
41   PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE));
42   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
43 
44   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
45 
46   PetscCall(DMGetLocalSection(*dm_sgs, &section));
47   PetscCall(PetscSectionSetFieldName(section, 0, ""));
48   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
49   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
50   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
51   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
52   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
53   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
54   PetscFunctionReturn(PETSC_SUCCESS);
55 };
56 
57 // @brief Evaluate data-driven SGS using fused method
58 static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
59   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
60   PetscMemType q_mem_type;
61 
62   PetscFunctionBeginUser;
63   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
64 
65   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
66 
67   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
68   PetscFunctionReturn(PETSC_SUCCESS);
69 }
70 
71 // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
72 static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
73   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
74   CeedQFunction       qf_sgs_dd_nodal;
75   CeedOperator        op_sgs_dd_nodal;
76   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
77   PetscInt            dim;
78   CeedVector          inv_multiplicity;
79   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
80   DMLabel             domain_label = NULL;
81   PetscInt            label_value = 0, height = 0, dm_field = 0;
82 
83   PetscFunctionBeginUser;
84   PetscCall(DMGetDimension(honee->dm, &dim));
85   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
86   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
87   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
88 
89   {  // Get velocity gradient information
90     CeedOperatorField op_field;
91     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
92     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
93     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
94   }
95   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
96   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
97 
98   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
99                                    &inv_multiplicity));
100 
101   // -- Create operator for SGS DD model nodal evaluation
102   switch (honee->phys->state_var) {
103     case STATEVAR_PRIMITIVE:
104       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
105       break;
106     case STATEVAR_CONSERVATIVE:
107       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
108       break;
109     case STATEVAR_ENTROPY:
110       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
111       break;
112   }
113 
114   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
115   CeedBasis basis_x_to_q;
116   PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q));
117 
118   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
119   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
120   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
121   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
122   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
123   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
124   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
125 
126   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
127   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
128   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord));
129   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
130   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
131                                            sgs_dd_setup_data->grid_aniso_ceed));
132   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
133   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
134 
135   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,
136                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
137 
138   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
139   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
140 
141   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
142   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
143   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
144   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
145   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
146   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
147   PetscFunctionReturn(PETSC_SUCCESS);
148 }
149 
150 // @brief Setup data-driven model inference using libCEED native implementation
151 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
152                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
153                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
154                                                                 void **ctx) {
155   CeedQFunction         qf_sgs_dd_inference;
156   CeedOperator          op_sgs_dd_inference;
157   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
158 
159   PetscFunctionBeginUser;
160   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
161                                                   &qf_sgs_dd_inference));
162 
163   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
164   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
165   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
166   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
167 
168   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
169   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
170   PetscCallCeed(ceed,
171                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
172   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
173 
174   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
175                                        op_context));
176   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
177 
178   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
179   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
180   PetscFunctionReturn(PETSC_SUCCESS);
181 }
182 
183 // @brief Perform data-driven model inference using libCEED native implementation
184 PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
185   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
186 
187   PetscFunctionBeginUser;
188   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
189   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
190   PetscCall(PetscLogGpuTimeBegin());
191   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
192   PetscCall(PetscLogGpuTimeEnd());
193   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
194   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
195   PetscFunctionReturn(PETSC_SUCCESS);
196 }
197 
198 // @brief Setup data-driven model inference using libtorch
199 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
200                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
201                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
202                                                                  void **ctx) {
203   const char     *ceed_resource;
204   char            model_path[PETSC_MAX_PATH_LEN] = "";
205   TorchDeviceType model_device_type;
206 
207   PetscFunctionBeginUser;
208   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
209   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
210   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
211   // On-device XPU is not working reliably currently, default to CPU inference evaluation
212   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
213   else model_device_type = TORCH_DEVICE_CPU;
214   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
215   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
216 
217   PetscCall(LoadModel_Torch(model_path, model_device_type));
218 
219   PetscFunctionReturn(PETSC_SUCCESS);
220 }
221 
222 // @brief Perform data-driven model inference using libtorch
223 static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
224   static PetscBool run_through = PETSC_FALSE;
225   PetscFunctionBeginUser;
226   if (!run_through) {
227     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
228   }
229   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
230   if (!run_through) {
231     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
232     run_through = PETSC_TRUE;
233   }
234   PetscFunctionReturn(PETSC_SUCCESS);
235 }
236 
237 // @brief Evaluate data-driven SGS using sequential method
238 PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
239   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
240   PetscMemType q_mem_type;
241   Vec          DD_Inputs_loc, DD_Outputs_loc;
242 
243   PetscFunctionBeginUser;
244   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
245   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
246   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
247 
248   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
249   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
250   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
251 
252   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
253   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
254   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
255   PetscFunctionReturn(PETSC_SUCCESS);
256 }
257 
258 // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
259 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
260   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
261   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
262   PetscInt            dim;
263   CeedVector          inv_multiplicity, eigvec;
264   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
265       elem_restr_dd_outputs;
266   DMLabel  domain_label = NULL;
267   PetscInt label_value = 0, height = 0, dm_field = 0;
268 
269   PetscFunctionBeginUser;
270   {  // Create DMs for data-driven input and output values
271     PetscSection section;
272     PetscInt     degree, q_extra;
273     {  // Get degree and number of quadrature points from dm_sgs
274       PetscFE         fe;
275       PetscSpace      basis;
276       PetscQuadrature quadrature;
277       PetscInt        num_qpnts;
278       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
279       PetscCall(PetscFEGetBasisSpace(fe, &basis));
280       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
281       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
282       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
283       q_extra = degree - num_qpnts;
284     }
285 
286     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
287     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE));
288     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
289     PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_inputs, sgs_dd_data->dm_dd_inputs));
290     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
291     PetscCall(PetscSectionSetFieldName(section, 0, ""));
292     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
293       char component_name[PETSC_MAX_PATH_LEN];
294 
295       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
296       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
297     }
298 
299     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
300     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE));
301     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
302     PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_outputs, sgs_dd_data->dm_dd_outputs));
303     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
304     PetscCall(PetscSectionSetFieldName(section, 0, ""));
305     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
306       char component_name[PETSC_MAX_PATH_LEN];
307 
308       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
309       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
310     }
311   }
312 
313   PetscCall(DMGetDimension(honee->dm, &dim));
314   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
315   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
316   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
317 
318   {  // Get velocity gradient information
319     CeedOperatorField op_field;
320     PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
321     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
322     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
323     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
324   }
325 
326   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
327   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
328   PetscCall(
329       DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec));
330   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
331 
332   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs));
333   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs));
334 
335   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
336                                    &inv_multiplicity));
337 
338   {  // Create operator for data-driven input evaluation
339     CeedQFunction qf_sgs_dd_inputs;
340     CeedOperator  op_sgs_dd_inputs;
341 
342     switch (honee->phys->state_var) {
343       case STATEVAR_PRIMITIVE:
344         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
345                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
346         break;
347       case STATEVAR_CONSERVATIVE:
348         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
349                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
350         break;
351       case STATEVAR_ENTROPY:
352         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
353                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
354         break;
355     }
356 
357     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
358     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
359     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
360     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
361     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
362     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
363     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
364 
365     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
366     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
367     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
368     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
369                                              sgs_dd_setup_data->grid_aniso_ceed));
370     PetscCallCeed(ceed,
371                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
372     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
373     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
374 
375     PetscCall(OperatorApplyContextCreate(honee->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
376                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
377     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
378     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
379   }
380 
381   {  // Create operator for data-driven output handling
382     CeedQFunction qf_sgs_dd_outputs;
383     CeedOperator  op_sgs_dd_outputs;
384 
385     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
386                                                     &qf_sgs_dd_outputs));
387     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
388     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
389     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
390     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
391     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
392     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
393 
394     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
395     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
396     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
397                                              sgs_dd_setup_data->grid_aniso_ceed));
398     PetscCallCeed(ceed,
399                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
400     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
401     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
402 
403     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,
404                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
405     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
406     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
407   }
408 
409   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
410 
411   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
412     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
413     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
414                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
415   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
416     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
417     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
418                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
419   }
420 
421   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
422 
423   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
424   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
425   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
426   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
427   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
428   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
429   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
430   PetscFunctionReturn(PETSC_SUCCESS);
431 }
432 
433 // @brief Create CeedOperator to compute SGS contribution to the residual
434 static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
435   SgsDDData           sgs_dd_data = honee->sgs_dd_data;
436   CeedInt             num_comp_q, q_data_size, num_comp_x;
437   PetscInt            dim;
438   CeedQFunction       qf_sgs_apply;
439   CeedOperator        op_sgs_apply;
440   CeedBasis           basis_sgs;
441   CeedVector          q_data;
442   CeedElemRestriction elem_restr_qd;
443 
444   PetscFunctionBeginUser;
445   PetscCall(DMGetDimension(honee->dm, &dim));
446   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q));
447   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
448 
449   {
450     DMLabel  domain_label = NULL;
451     PetscInt label_value = 0, height = 0, dm_field = 0;
452 
453     PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs));
454     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,
455                        &q_data, &q_data_size));
456   }
457 
458   switch (honee->phys->state_var) {
459     case STATEVAR_PRIMITIVE:
460       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
461       break;
462     case STATEVAR_CONSERVATIVE:
463       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
464       break;
465     case STATEVAR_ENTROPY:
466       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
467       break;
468   }
469 
470   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
471   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
472   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE));
473   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
474   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
475 
476   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
477   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
478   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
479   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
480   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
481 
482   PetscCall(
483       OperatorApplyContextCreate(honee->dm, honee->dm, ceed, op_sgs_apply, honee->q_ceed, honee->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx));
484 
485   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
486   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
487   PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs));
488   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
489   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
490   PetscFunctionReturn(PETSC_SUCCESS);
491 }
492 
493 // @brief Calculate and add data-driven SGS residual to the global residual
494 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) {
495   SgsDDData    sgs_dd_data = honee->sgs_dd_data;
496   Vec          VelocityGradient, SGSNodal_loc;
497   PetscMemType sgs_nodal_mem_type;
498 
499   PetscFunctionBeginUser;
500   PetscCall(PetscLogEventBegin(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
501   PetscCall(DMGetGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
502   PetscCall(VelocityGradientProjectionApply(honee->grad_velo_proj, Q_loc, VelocityGradient));
503 
504   // -- Compute Nodal SGS tensor
505   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
506   PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc));
507 
508   // -- Compute contribution of the SGS stress
509   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
510   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
511 
512   // -- Return local SGS vector
513   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
514   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
515   PetscCall(DMRestoreGlobalVector(honee->grad_velo_proj->dm, &VelocityGradient));
516   PetscCall(PetscLogEventEnd(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
517   PetscFunctionReturn(PETSC_SUCCESS);
518 }
519 
520 // @brief B = A^T, A is NxM, B is MxN
521 static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
522   PetscFunctionBeginUser;
523   for (PetscInt i = 0; i < N; i++) {
524     for (PetscInt j = 0; j < M; j++) {
525       B[j * N + i] = A[i * M + j];
526     }
527   }
528   PetscFunctionReturn(PETSC_SUCCESS);
529 }
530 
531 // @brief Read neural network coefficients from file and put into context struct
532 static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
533   SgsDDContext sgsdd_ctx;
534   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
535   char         file_path[PETSC_MAX_PATH_LEN];
536   PetscScalar *temp;
537 
538   PetscFunctionBeginUser;
539   {
540     SgsDDContext sgsdd_temp;
541     PetscCall(PetscNew(&sgsdd_temp));
542     *sgsdd_temp                     = **psgsdd_ctx;
543     sgsdd_temp->offsets.bias1       = 0;
544     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
545     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
546     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
547     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
548     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
549     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
550     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
551     *sgsdd_ctx = *sgsdd_temp;
552     PetscCall(PetscFree(sgsdd_temp));
553   }
554 
555   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
556   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
557   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
558   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
559   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
560   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
561 
562   {
563     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
564     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
565     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
566     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
567     PetscCall(PetscFree(temp));
568   }
569   {
570     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
571     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
572     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
573     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
574     PetscCall(PetscFree(temp));
575   }
576 
577   PetscCall(PetscFree(*psgsdd_ctx));
578   *psgsdd_ctx = sgsdd_ctx;
579   PetscFunctionReturn(PETSC_SUCCESS);
580 }
581 
582 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) {
583   PetscReal                alpha = 0;
584   SgsDDContext             sgsdd_ctx;
585   MPI_Comm                 comm                           = honee->comm;
586   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
587   SgsDDSetupData           sgs_dd_setup_data;
588   NewtonianIdealGasContext gas;
589 
590   PetscFunctionBeginUser;
591   PetscCall(
592       VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &honee->grad_velo_proj));
593 
594   PetscCall(PetscNew(&honee->sgs_dd_data));
595   honee->sgs_dd_data->num_comp_inputs  = 6;
596   honee->sgs_dd_data->num_comp_outputs = 6;
597 
598   PetscCall(PetscNew(&sgs_dd_setup_data));
599 
600   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
601   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
602   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
603                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
604   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
605   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
606   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
607                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
608                              NULL));
609   PetscOptionsEnd();
610 
611   PetscCall(PetscNew(&sgsdd_ctx));
612   sgsdd_ctx->num_layers  = 1;
613   sgsdd_ctx->num_inputs  = 6;
614   sgsdd_ctx->num_outputs = 6;
615   sgsdd_ctx->num_neurons = 20;
616   sgsdd_ctx->alpha       = alpha;
617 
618   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
619 
620   // -- Create DM for storing SGS tensor at nodes
621   PetscCall(
622       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));
623 
624   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
625   sgsdd_ctx->gas = *gas;
626   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
627   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
628   PetscCallCeed(ceed,
629                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
630   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
631 
632   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx));
633 
634   // -- Compute and store anisotropy tensor
635   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed));
636 
637   // -- Create Nodal Evaluation Operator
638   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
639     case SGS_MODEL_DD_FUSED:
640       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data));
641       break;
642     case SGS_MODEL_DD_SEQENTIAL_CEED:
643     case SGS_MODEL_DD_SEQENTIAL_TORCH:
644       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data));
645       break;
646   }
647 
648   // -- Create Operator to evalutate residual of SGS stress
649   PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data));
650 
651   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
652   PetscFunctionReturn(PETSC_SUCCESS);
653 }
654 
655 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
656   PetscFunctionBeginUser;
657   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
658   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
659 
660   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
661   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed));
662   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
663   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
664   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx));
665   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx));
666   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
667   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs));
668   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs));
669   if (sgs_dd_data->sgs_nodal_inference_ctx) PetscCall(sgs_dd_data->sgs_nodal_inference_ctx_destroy(sgs_dd_data->sgs_nodal_inference_ctx));
670   PetscCall(PetscFree(sgs_dd_data));
671   PetscFunctionReturn(PETSC_SUCCESS);
672 }
673