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