xref: /honee/problems/sgs_dd_model.c (revision ad494f68b826069604e1f8028d8ad837eec0f9a5)
1 // Copyright (c) 2017-2023, 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 "../navierstokes.h"
13 
14 typedef struct {
15   CeedElemRestriction  elem_restr_grid_aniso, elem_restr_sgs;
16   CeedVector           grid_aniso_ceed;
17   CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx;
18 } *SgsDDSetupData;
19 
20 PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
21   Ceed ceed;
22 
23   PetscFunctionBeginUser;
24   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
25 
26   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
27   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
28   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
29   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
30   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
31   PetscCall(PetscFree(sgs_dd_setup_data));
32   PetscFunctionReturn(PETSC_SUCCESS);
33 }
34 
35 // @brief Create DM for storing subgrid stress at nodes
36 static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
37   PetscSection section;
38 
39   PetscFunctionBeginUser;
40   *num_components = 6;
41 
42   PetscCall(DMClone(dm_source, dm_sgs));
43   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
44 
45   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
46 
47   PetscCall(DMGetLocalSection(*dm_sgs, &section));
48   PetscCall(PetscSectionSetFieldName(section, 0, ""));
49   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
50   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
51   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
52   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
53   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
54   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
55   PetscFunctionReturn(PETSC_SUCCESS);
56 };
57 
58 // @brief Evaluate data-driven SGS using fused method
59 static PetscErrorCode SgsDDNodalStressEval_Fused(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
60   SgsDDData    sgs_dd_data = user->sgs_dd_data;
61   PetscMemType q_mem_type;
62 
63   PetscFunctionBeginUser;
64   PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
65 
66   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
67 
68   PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc));
69   PetscFunctionReturn(PETSC_SUCCESS);
70 }
71 
72 // @brief Create CeedOperator to calculate data-drive SGS at nodes
73 static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
74   SgsDDData           sgs_dd_data = user->sgs_dd_data;
75   CeedQFunction       qf_multiplicity, qf_sgs_dd_nodal;
76   CeedOperator        op_multiplicity, op_sgs_dd_nodal;
77   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
78   PetscInt            dim;
79   CeedVector          multiplicity, inv_multiplicity;
80   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
81   DMLabel             domain_label = NULL;
82   PetscInt            label_value = 0, height = 0, dm_field = 0;
83 
84   PetscFunctionBeginUser;
85   PetscCall(DMGetDimension(user->dm, &dim));
86   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
87   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
88   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
89 
90   {  // Get velocity gradient information
91     CeedOperatorField op_field;
92     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
93     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
94     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
95   }
96   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
97   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
98 
99   // -- Create inverse multiplicity for correcting nodal assembly
100   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL));
101   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity));
102   PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, 1, &elem_restr_inv_multiplicity));
103   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL));
104 
105   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity));
106   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE));
107   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE));
108 
109   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity));
110   PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling"));
111   PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
112   PetscCallCeed(ceed,
113                 CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
114 
115   PetscCallCeed(ceed, CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE));
116 
117   // -- Create operator for SGS DD model nodal evaluation
118   switch (user->phys->state_var) {
119     case STATEVAR_PRIMITIVE:
120       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
121       break;
122     case STATEVAR_CONSERVATIVE:
123       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
124       break;
125     default:
126       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Data-driven SGS nodal evaluation not available for chosen state variable");
127   }
128 
129   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
130   CeedBasis basis_x_to_q;
131   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
132 
133   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
134   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
135   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
136   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
137   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
138   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
139   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
140 
141   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
142   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
143   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord));
144   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
145   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
146                                            sgs_dd_setup_data->grid_aniso_ceed));
147   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
148   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
149 
150   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,
151                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
152 
153   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
154   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
155 
156   PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity));
157   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
158   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
159   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
160   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity));
161   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
162   PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity));
163   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
164   PetscFunctionReturn(PETSC_SUCCESS);
165 }
166 
167 // @brief Create CeedOperator to compute SGS contribution to the residual
168 static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
169   SgsDDData     sgs_dd_data = user->sgs_dd_data;
170   CeedInt       num_comp_q, num_comp_qd, num_comp_x;
171   PetscInt      dim;
172   CeedQFunction qf_sgs_apply;
173   CeedOperator  op_sgs_apply;
174   CeedBasis     basis_sgs;
175 
176   PetscFunctionBeginUser;
177   PetscCall(DMGetDimension(user->dm, &dim));
178   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
179   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd));
180   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
181 
182   PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs));
183 
184   switch (user->phys->state_var) {
185     case STATEVAR_PRIMITIVE:
186       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
187       break;
188     case STATEVAR_CONSERVATIVE:
189       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
190       break;
191     default:
192       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Nodal SGS evaluation not available for chosen state variable");
193   }
194 
195   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
196   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
197   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE));
198   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
199   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
200 
201   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
202   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
203   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
204   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
205   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
206 
207   PetscCall(
208       OperatorApplyContextCreate(user->dm, user->dm, ceed, op_sgs_apply, user->q_ceed, user->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx));
209 
210   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
211   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
212   PetscFunctionReturn(PETSC_SUCCESS);
213 }
214 
215 // @brief Calculate and add data-driven SGS residual to the global residual
216 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc) {
217   SgsDDData    sgs_dd_data = user->sgs_dd_data;
218   Vec          VelocityGradient, SGSNodal_loc;
219   PetscMemType sgs_nodal_mem_type;
220 
221   PetscFunctionBeginUser;
222   PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
223   PetscCall(VelocityGradientProjectionApply(user->grad_velo_proj, Q_loc, VelocityGradient));
224 
225   // -- Compute Nodal SGS tensor
226   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
227   PetscCall(sgs_dd_data->sgs_nodal_eval(user, Q_loc, VelocityGradient, SGSNodal_loc));
228 
229   // -- Compute contribution of the SGS stress
230   PetscCall(VecP2C(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
231   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
232 
233   // -- Return local SGS vector
234   PetscCall(VecC2P(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
235   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
236   PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
237   PetscFunctionReturn(PETSC_SUCCESS);
238 }
239 
240 // @brief B = A^T, A is NxM, B is MxN
241 static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
242   PetscFunctionBeginUser;
243   for (PetscInt i = 0; i < N; i++) {
244     for (PetscInt j = 0; j < M; j++) {
245       B[j * N + i] = A[i * M + j];
246     }
247   }
248   PetscFunctionReturn(PETSC_SUCCESS);
249 }
250 
251 // @brief Read neural network coefficients from file and put into context struct
252 static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
253   SgsDDContext sgsdd_ctx;
254   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
255   char         file_path[PETSC_MAX_PATH_LEN];
256   PetscScalar *temp;
257 
258   PetscFunctionBeginUser;
259   {
260     SgsDDContext sgsdd_temp;
261     PetscCall(PetscNew(&sgsdd_temp));
262     *sgsdd_temp                     = **psgsdd_ctx;
263     sgsdd_temp->offsets.bias1       = 0;
264     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
265     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
266     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
267     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
268     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
269     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
270     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
271     *sgsdd_ctx = *sgsdd_temp;
272     PetscCall(PetscFree(sgsdd_temp));
273   }
274 
275   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
276   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
277   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
278   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
279   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
280   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
281 
282   {
283     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
284     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
285     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
286     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
287     PetscCall(PetscFree(temp));
288   }
289   {
290     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
291     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
292     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
293     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
294     PetscCall(PetscFree(temp));
295   }
296 
297   PetscCall(PetscFree(*psgsdd_ctx));
298   *psgsdd_ctx = sgsdd_ctx;
299   PetscFunctionReturn(PETSC_SUCCESS);
300 }
301 
302 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
303   PetscReal                alpha = 0;
304   SgsDDContext             sgsdd_ctx;
305   MPI_Comm                 comm                           = user->comm;
306   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
307   SgsDDSetupData           sgs_dd_setup_data;
308   NewtonianIdealGasContext gas;
309 
310   PetscFunctionBeginUser;
311   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, user->phys->state_var, ceed_data->elem_restr_q, ceed_data->basis_q,
312                                             &user->grad_velo_proj));
313 
314   PetscCall(PetscNew(&sgsdd_ctx));
315 
316   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
317   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
318   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
319                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
320   PetscOptionsEnd();
321 
322   sgsdd_ctx->num_layers  = 1;
323   sgsdd_ctx->num_inputs  = 6;
324   sgsdd_ctx->num_outputs = 6;
325   sgsdd_ctx->num_neurons = 20;
326   sgsdd_ctx->alpha       = alpha;
327 
328   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
329 
330   // -- Create DM for storing SGS tensor at nodes
331   PetscCall(PetscNew(&user->sgs_dd_data));
332   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));
333 
334   PetscCall(PetscNew(&sgs_dd_setup_data));
335 
336   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
337   sgsdd_ctx->gas = *gas;
338   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
339   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
340   PetscCallCeed(ceed,
341                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
342   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
343 
344   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfunction_context, &sgs_dd_setup_data->ifunction_qfctx));
345 
346   // -- Compute and store anisotropy tensor
347   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
348                                                      &sgs_dd_setup_data->grid_aniso_ceed));
349 
350   // -- Create Nodal Evaluation Operator
351   PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, user, ceed_data, sgs_dd_setup_data));
352 
353   // -- Create Operator to evalutate residual of SGS stress
354   PetscCall(SgsSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data));
355 
356   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
357   PetscFunctionReturn(PETSC_SUCCESS);
358 }
359 
360 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
361   PetscFunctionBeginUser;
362   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
363   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
364 
365   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
366   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
367   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
368   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
369   PetscCall(PetscFree(sgs_dd_data));
370   PetscFunctionReturn(PETSC_SUCCESS);
371 }
372