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