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