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