xref: /honee/problems/sgs_dd_model.c (revision b4c37c5c26bd208d7f474cbd9aa0850e82e1a854)
162b7942eSJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors.
262b7942eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
362b7942eSJames Wright //
462b7942eSJames Wright // SPDX-License-Identifier: BSD-2-Clause
562b7942eSJames Wright //
662b7942eSJames Wright // This file is part of CEED:  http://github.com/ceed
762b7942eSJames Wright 
862b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h"
962b7942eSJames Wright 
1062b7942eSJames Wright #include <petscdmplex.h>
1162b7942eSJames Wright 
1262b7942eSJames Wright #include "../navierstokes.h"
1362b7942eSJames Wright 
14c38c977aSJames Wright typedef struct {
15c38c977aSJames Wright   CeedElemRestriction  elem_restr_grid_aniso, elem_restr_sgs;
16c38c977aSJames Wright   CeedVector           grid_aniso_ceed;
17ee1455b7SJames Wright   CeedQFunctionContext sgsdd_qfctx;
18c38c977aSJames Wright } *SGS_DD_ModelSetupData;
19c38c977aSJames Wright 
20c38c977aSJames Wright PetscErrorCode SGS_DD_ModelSetupDataDestroy(SGS_DD_ModelSetupData sgs_dd_setup_data) {
21*b4c37c5cSJames Wright   Ceed ceed;
22*b4c37c5cSJames Wright 
23c38c977aSJames Wright   PetscFunctionBeginUser;
24*b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
25*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
26*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
27*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
28*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
29c38c977aSJames Wright 
30c38c977aSJames Wright   PetscCall(PetscFree(sgs_dd_setup_data));
31d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32c38c977aSJames Wright }
33c38c977aSJames Wright 
34ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
3567263decSKenneth E. Jansen PetscErrorCode SGS_DD_ModelCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
36ee1455b7SJames Wright   PetscFE      fe;
37ee1455b7SJames Wright   PetscSection section;
38ee1455b7SJames Wright   PetscInt     dim;
39ee1455b7SJames Wright 
40ee1455b7SJames Wright   PetscFunctionBeginUser;
41ee1455b7SJames Wright   *num_components  = 6;
4267263decSKenneth E. Jansen   PetscInt q_order = degree + q_extra;
43ee1455b7SJames Wright 
44ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
45ee1455b7SJames Wright   PetscCall(DMGetDimension(*dm_sgs, &dim));
46ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
47ee1455b7SJames Wright 
4867263decSKenneth E. Jansen   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, *num_components, PETSC_FALSE, degree, q_order, &fe));
49ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "Subgrid Stress Projection"));
50ee1455b7SJames Wright   PetscCall(DMAddField(*dm_sgs, NULL, (PetscObject)fe));
51ee1455b7SJames Wright   PetscCall(DMCreateDS(*dm_sgs));
52ee1455b7SJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(*dm_sgs, PETSC_DETERMINE, NULL));
53ee1455b7SJames Wright 
54ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
55ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
56ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
57ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
58ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
59ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
60ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
61ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
62ee1455b7SJames Wright 
63ee1455b7SJames Wright   PetscCall(PetscFEDestroy(&fe));
64ee1455b7SJames Wright 
65d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66ee1455b7SJames Wright };
67ee1455b7SJames Wright 
68ee1455b7SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes
69ee1455b7SJames Wright PetscErrorCode SGS_DD_ModelSetupNodalEvaluation(Ceed ceed, User user, CeedData ceed_data, SGS_DD_ModelSetupData sgs_dd_setup_data) {
70ee1455b7SJames Wright   SGS_DD_Data         sgs_dd_data = user->sgs_dd_data;
71ee1455b7SJames Wright   CeedQFunction       qf_multiplicity, qf_sgs_dd_nodal;
72ee1455b7SJames Wright   CeedOperator        op_multiplicity, op_sgs_dd_nodal;
7367263decSKenneth E. Jansen   CeedInt             num_elem, elem_size, num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
74defe8520SJames Wright   PetscInt            dim;
754b0f6111SJames Wright   CeedVector          multiplicity, inv_multiplicity;
76ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
77ee1455b7SJames Wright 
78ee1455b7SJames Wright   PetscFunctionBeginUser;
79ee1455b7SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
80*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
81*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
82*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
83*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &num_elem));
84*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetElementSize(ceed_data->elem_restr_q, &elem_size));
85ee1455b7SJames Wright 
86ee1455b7SJames Wright   {  // Get velocity gradient information
87ee1455b7SJames Wright     CeedOperatorField op_field;
88*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
89*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
90*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
91ee1455b7SJames Wright   }
9267263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, -1, 0, &elem_restr_sgs, NULL, NULL));
93*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
94ee1455b7SJames Wright 
95ee1455b7SJames Wright   // -- Create inverse multiplicity for correcting nodal assembly
96*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL));
97*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity));
98*b4c37c5cSJames Wright   PetscCallCeed(
99*b4c37c5cSJames Wright       ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, 1, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_inv_multiplicity));
100*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL));
101ee1455b7SJames Wright 
102*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity));
103*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE));
104*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE));
105ee1455b7SJames Wright 
106*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity));
107*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling"));
108*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
109*b4c37c5cSJames Wright   PetscCallCeed(
110*b4c37c5cSJames Wright       ceed, CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
111ee1455b7SJames Wright 
112*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE));
113ee1455b7SJames Wright 
114ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
115ee1455b7SJames Wright   switch (user->phys->state_var) {
116ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
117*b4c37c5cSJames Wright       PetscCallCeed(
118*b4c37c5cSJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Prim, ComputeSGS_DDAnisotropicNodal_Prim_loc, &qf_sgs_dd_nodal));
119ee1455b7SJames Wright       break;
120ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
121*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicNodal_Conserv, ComputeSGS_DDAnisotropicNodal_Conserv_loc,
122*b4c37c5cSJames Wright                                                       &qf_sgs_dd_nodal));
123ee1455b7SJames Wright       break;
124ee1455b7SJames Wright     default:
125ee1455b7SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP,
126ee1455b7SJames Wright               "Anisotropic data-driven SGS nodal evaluation not available for chosen state variable");
127ee1455b7SJames Wright   }
128ee1455b7SJames Wright 
129ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
130ee1455b7SJames Wright   CeedBasis basis_x_to_q;
131*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
132ee1455b7SJames Wright 
133*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
134*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
135*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
136*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
137*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
138*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
139*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
140ee1455b7SJames Wright 
141*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
142*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, user->q_ceed));
143*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord));
144*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
145*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_COLLOCATED,
146*b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
147*b4c37c5cSJames Wright   PetscCallCeed(ceed,
148*b4c37c5cSJames Wright                 CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, inv_multiplicity));
149*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
150ee1455b7SJames Wright 
1514b0f6111SJames Wright   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,
1524b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
153ee1455b7SJames Wright 
154ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
155ee1455b7SJames Wright 
156*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity));
157*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
158*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
159*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
160*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity));
161*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
162*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity));
163*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
164d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
165ee1455b7SJames Wright }
166ee1455b7SJames Wright 
1679c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
1689c678832SJames Wright PetscErrorCode SGS_ModelSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SGS_DD_ModelSetupData sgs_dd_setup_data) {
1699c678832SJames Wright   SGS_DD_Data   sgs_dd_data = user->sgs_dd_data;
17067263decSKenneth E. Jansen   CeedInt       num_comp_q, num_comp_qd, num_comp_x;
171defe8520SJames Wright   PetscInt      dim;
1729c678832SJames Wright   CeedQFunction qf_sgs_apply;
1739c678832SJames Wright   CeedOperator  op_sgs_apply;
1749c678832SJames Wright   CeedBasis     basis_sgs;
1759c678832SJames Wright 
1769c678832SJames Wright   PetscFunctionBeginUser;
1779c678832SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
178*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
179*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd));
180*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
1819c678832SJames Wright 
18267263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs));
1839c678832SJames Wright 
1849c678832SJames Wright   switch (user->phys->state_var) {
1859c678832SJames Wright     case STATEVAR_PRIMITIVE:
186*b4c37c5cSJames Wright       PetscCallCeed(ceed,
187*b4c37c5cSJames Wright                     CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSubgridStress_Prim, IFunction_NodalSubgridStress_Prim_loc, &qf_sgs_apply));
1889c678832SJames Wright       break;
1899c678832SJames Wright     case STATEVAR_CONSERVATIVE:
190*b4c37c5cSJames Wright       PetscCallCeed(
191*b4c37c5cSJames Wright           ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSubgridStress_Conserv, IFunction_NodalSubgridStress_Conserv_loc, &qf_sgs_apply));
1929c678832SJames Wright       break;
1939c678832SJames Wright     default:
1949c678832SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Nodal SGS evaluation not available for chosen state variable");
1959c678832SJames Wright   }
1969c678832SJames Wright 
197*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->sgsdd_qfctx));
198*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
199*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE));
200*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "x", num_comp_x, CEED_EVAL_INTERP));
201*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
202*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
2039c678832SJames Wright 
204*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
205*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
206*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
207*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
208*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
209*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
2109c678832SJames Wright 
2119c678832SJames Wright   PetscCall(
2129c678832SJames Wright       OperatorApplyContextCreate(user->dm, user->dm, ceed, op_sgs_apply, user->q_ceed, user->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx));
2139c678832SJames Wright 
214*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
215*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
216d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2179c678832SJames Wright }
2189c678832SJames Wright 
2199c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
2209c678832SJames Wright PetscErrorCode SGS_DD_ModelApplyIFunction(User user, const Vec Q_loc, Vec G_loc) {
2219c678832SJames Wright   SGS_DD_Data  sgs_dd_data = user->sgs_dd_data;
2229c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
2239c678832SJames Wright   PetscMemType sgs_nodal_mem_type, q_mem_type;
2249c678832SJames Wright 
2259c678832SJames Wright   PetscFunctionBeginUser;
2269c678832SJames Wright   PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
2279c678832SJames Wright   PetscCall(VelocityGradientProjectionApply(user, Q_loc, VelocityGradient));
2289c678832SJames Wright 
2299c678832SJames Wright   // -- Compute Nodal SGS tensor
2309c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
2319c678832SJames Wright   PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
2329c678832SJames Wright 
2339c678832SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
2349c678832SJames Wright 
2359c678832SJames Wright   PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc));
2369c678832SJames Wright   PetscCall(VecP2C(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
2379c678832SJames Wright 
2389c678832SJames Wright   // -- Compute contribution of the SGS stress
2399c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
2409c678832SJames Wright 
2419c678832SJames Wright   // -- Return local SGS vector
2429c678832SJames Wright   PetscCall(VecC2P(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
2439c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
2449c678832SJames Wright   PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
2459c678832SJames Wright 
246d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2479c678832SJames Wright }
2489c678832SJames Wright 
24962b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
25062b7942eSJames Wright PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
25162b7942eSJames Wright   PetscFunctionBeginUser;
25262b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
25362b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
25462b7942eSJames Wright       B[j * N + i] = A[i * M + j];
25562b7942eSJames Wright     }
25662b7942eSJames Wright   }
257d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
25862b7942eSJames Wright }
25962b7942eSJames Wright 
26062b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
26162b7942eSJames Wright PetscErrorCode SGS_DD_ModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SGS_DDModelContext *psgsdd_ctx) {
26262b7942eSJames Wright   SGS_DDModelContext sgsdd_ctx;
26362b7942eSJames Wright   PetscInt           num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
26462b7942eSJames Wright   char               file_path[PETSC_MAX_PATH_LEN];
26562b7942eSJames Wright   PetscScalar       *temp;
26662b7942eSJames Wright 
26762b7942eSJames Wright   PetscFunctionBeginUser;
26862b7942eSJames Wright   {
26962b7942eSJames Wright     SGS_DDModelContext sgsdd_temp;
27062b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
27162b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
27262b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
27362b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
27462b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
27562b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
27662b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
27762b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
27862b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
27962b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
28062b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
28162b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
28262b7942eSJames Wright   }
28362b7942eSJames Wright 
28462b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
28562b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
28662b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
28762b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
28862b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
28962b7942eSJames Wright   PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
29062b7942eSJames Wright 
29162b7942eSJames Wright   {
29262b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
29362b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
29462b7942eSJames Wright     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
29562b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
29662b7942eSJames Wright     PetscCall(PetscFree(temp));
29762b7942eSJames Wright   }
29862b7942eSJames Wright   {
29962b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
30062b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
30162b7942eSJames Wright     PetscCall(PHASTADatFileReadToArrayReal(comm, file_path, temp));
30262b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
30362b7942eSJames Wright     PetscCall(PetscFree(temp));
30462b7942eSJames Wright   }
30562b7942eSJames Wright 
30662b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
30762b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
308d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
30962b7942eSJames Wright }
31062b7942eSJames Wright 
31162b7942eSJames Wright PetscErrorCode SGS_DD_ModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
312ee1455b7SJames Wright   PetscReal                alpha = 0;
31362b7942eSJames Wright   SGS_DDModelContext       sgsdd_ctx;
31462b7942eSJames Wright   MPI_Comm                 comm                           = user->comm;
315ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
316c38c977aSJames Wright   SGS_DD_ModelSetupData    sgs_dd_setup_data;
317ee1455b7SJames Wright   NewtonianIdealGasContext gas;
31862b7942eSJames Wright   PetscFunctionBeginUser;
31962b7942eSJames Wright 
3209ab09d51SJames Wright   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem));
3219ab09d51SJames Wright 
32262b7942eSJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
32362b7942eSJames Wright 
3244b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
32562b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
32662b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
32762b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
32862b7942eSJames Wright   PetscOptionsEnd();
32962b7942eSJames Wright 
3304b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
33162b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
33262b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
33362b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
33462b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
33562b7942eSJames Wright 
33601ab89c1SJames Wright   PetscCall(SGS_DD_ModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
33762b7942eSJames Wright 
338ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
339ee1455b7SJames Wright   PetscCall(PetscNew(&user->sgs_dd_data));
34067263decSKenneth E. Jansen   PetscCall(
34167263decSKenneth E. Jansen       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));
342ee1455b7SJames Wright 
343c38c977aSJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
344c38c977aSJames Wright 
345*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
346ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
347*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
348*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
349*b4c37c5cSJames Wright   PetscCallCeed(ceed,
350*b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
351*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
352ee1455b7SJames Wright 
353c38c977aSJames Wright   // -- Compute and store anisotropy tensor
354c38c977aSJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
355c38c977aSJames Wright                                                      &sgs_dd_setup_data->grid_aniso_ceed));
356c38c977aSJames Wright 
357ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
358ee1455b7SJames Wright   PetscCall(SGS_DD_ModelSetupNodalEvaluation(ceed, user, ceed_data, sgs_dd_setup_data));
359ee1455b7SJames Wright 
3609c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
3619c678832SJames Wright   PetscCall(SGS_ModelSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data));
3629c678832SJames Wright 
363c38c977aSJames Wright   PetscCall(SGS_DD_ModelSetupDataDestroy(sgs_dd_setup_data));
364d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
36562b7942eSJames Wright }
366ee1455b7SJames Wright 
367ee1455b7SJames Wright PetscErrorCode SGS_DD_DataDestroy(SGS_DD_Data sgs_dd_data) {
368ee1455b7SJames Wright   PetscFunctionBeginUser;
369d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
370*b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
371ee1455b7SJames Wright 
372*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
373ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
374ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
375ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
376ee1455b7SJames Wright 
377d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
378ee1455b7SJames Wright }
379