xref: /honee/problems/sgs_dd_model.c (revision 4081638529d72d3eb1696fa2e9d59287ec507b73)
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;
17*40816385SJames Wright   CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx;
1842454adaSJames Wright } *SgsDDModelSetupData;
19c38c977aSJames Wright 
2042454adaSJames Wright PetscErrorCode SgsDDModelSetupDataDestroy(SgsDDModelSetupData sgs_dd_setup_data) {
21b4c37c5cSJames Wright   Ceed ceed;
22c38c977aSJames Wright   PetscFunctionBeginUser;
23b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
24b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
25b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
26b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
27b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
28*40816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
29c38c977aSJames Wright   PetscCall(PetscFree(sgs_dd_setup_data));
30d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
31c38c977aSJames Wright }
32c38c977aSJames Wright 
33ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
34cceb3143SJames Wright static PetscErrorCode SgsDDModelCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
35ee1455b7SJames Wright   PetscSection section;
36ee1455b7SJames Wright 
37ee1455b7SJames Wright   PetscFunctionBeginUser;
38ee1455b7SJames Wright   *num_components = 6;
39ee1455b7SJames Wright 
40ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
41ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
42ee1455b7SJames Wright 
43da4ca0cfSJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
44ee1455b7SJames Wright 
45ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
46ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
47ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
48ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
49ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
50ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
51ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
52ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
53d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
54ee1455b7SJames Wright };
55ee1455b7SJames Wright 
56cceb3143SJames Wright PetscErrorCode SgsDDModelNodalStressEval_Fused(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
57cceb3143SJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
58cceb3143SJames Wright   PetscMemType q_mem_type;
59cceb3143SJames Wright 
60cceb3143SJames Wright   PetscFunctionBeginUser;
61cceb3143SJames Wright   PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
62cceb3143SJames Wright 
63cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
64cceb3143SJames Wright 
65cceb3143SJames Wright   PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc));
66cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
67cceb3143SJames Wright }
68cceb3143SJames Wright 
69ee1455b7SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes
70cceb3143SJames Wright static PetscErrorCode SgsDDModelSetupNodalEvaluation_Fused(Ceed ceed, User user, CeedData ceed_data, SgsDDModelSetupData sgs_dd_setup_data) {
7142454adaSJames Wright   SgsDDData           sgs_dd_data = user->sgs_dd_data;
72ee1455b7SJames Wright   CeedQFunction       qf_multiplicity, qf_sgs_dd_nodal;
73ee1455b7SJames Wright   CeedOperator        op_multiplicity, op_sgs_dd_nodal;
7415c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
75defe8520SJames Wright   PetscInt            dim;
764b0f6111SJames Wright   CeedVector          multiplicity, inv_multiplicity;
77ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
7815c18037SJames Wright   DMLabel             domain_label = NULL;
7915c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
80ee1455b7SJames Wright 
81ee1455b7SJames Wright   PetscFunctionBeginUser;
82ee1455b7SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
83b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
84b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
85b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
86ee1455b7SJames Wright 
87ee1455b7SJames Wright   {  // Get velocity gradient information
88ee1455b7SJames Wright     CeedOperatorField op_field;
89b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
90b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
91b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
92ee1455b7SJames Wright   }
9315c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
94b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
95ee1455b7SJames Wright 
96ee1455b7SJames Wright   // -- Create inverse multiplicity for correcting nodal assembly
97b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL));
98b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity));
9915c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, 1, &elem_restr_inv_multiplicity));
100b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL));
101ee1455b7SJames Wright 
102b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity));
103b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE));
104b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE));
105ee1455b7SJames Wright 
106b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity));
107b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling"));
10858e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
10958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed,
11058e1cbfdSJeremy L Thompson                 CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
111ee1455b7SJames Wright 
112b4c37c5cSJames 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:
11742454adaSJames Wright       PetscCallCeed(ceed,
11842454adaSJames Wright                     CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDAnisotropicNodal_Prim, ComputeSgsDDAnisotropicNodal_Prim_loc, &qf_sgs_dd_nodal));
119ee1455b7SJames Wright       break;
120ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
12142454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDAnisotropicNodal_Conserv, ComputeSgsDDAnisotropicNodal_Conserv_loc,
122b4c37c5cSJames 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;
131b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
132ee1455b7SJames Wright 
133b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
134b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
135b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
136b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
137b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
138b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
139b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
140ee1455b7SJames Wright 
141b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
14258e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
143b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord));
14458e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
14558e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
146b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
14758e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
14858e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
149ee1455b7SJames Wright 
1504b0f6111SJames 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,
1514b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
152ee1455b7SJames Wright 
153ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
154cceb3143SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDModelNodalStressEval_Fused;
155ee1455b7SJames Wright 
156b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity));
157b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
158b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
159b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
160b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity));
161b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
162b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity));
163b4c37c5cSJames 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
168cceb3143SJames Wright static PetscErrorCode SgsModelSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDModelSetupData sgs_dd_setup_data) {
16942454adaSJames Wright   SgsDDData     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));
178b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
179b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd));
180b4c37c5cSJames 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:
18642454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
1879c678832SJames Wright       break;
1889c678832SJames Wright     case STATEVAR_CONSERVATIVE:
18942454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
1909c678832SJames Wright       break;
1919c678832SJames Wright     default:
1929c678832SJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Nodal SGS evaluation not available for chosen state variable");
1939c678832SJames Wright   }
1949c678832SJames Wright 
195*40816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
196b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
197b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE));
198b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
199b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
2009c678832SJames Wright 
201b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
202b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
20358e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
204b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
205b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
2069c678832SJames Wright 
2079c678832SJames Wright   PetscCall(
2089c678832SJames 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));
2099c678832SJames Wright 
210b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
211b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
212d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2139c678832SJames Wright }
2149c678832SJames Wright 
2159c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
21642454adaSJames Wright PetscErrorCode SgsDDModelApplyIFunction(User user, const Vec Q_loc, Vec G_loc) {
21742454adaSJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
2189c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
219cceb3143SJames Wright   PetscMemType sgs_nodal_mem_type;
2209c678832SJames Wright 
2219c678832SJames Wright   PetscFunctionBeginUser;
2229c678832SJames Wright   PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
223f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionApply(user->grad_velo_proj, Q_loc, VelocityGradient));
2249c678832SJames Wright 
2259c678832SJames Wright   // -- Compute Nodal SGS tensor
2269c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
227cceb3143SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(user, Q_loc, VelocityGradient, SGSNodal_loc));
2289c678832SJames Wright 
2299c678832SJames Wright   // -- Compute contribution of the SGS stress
230cceb3143SJames Wright   PetscCall(VecP2C(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
2319c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
2329c678832SJames Wright 
2339c678832SJames Wright   // -- Return local SGS vector
2349c678832SJames Wright   PetscCall(VecC2P(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
2359c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
2369c678832SJames Wright   PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
237d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2389c678832SJames Wright }
2399c678832SJames Wright 
24062b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
241cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
24262b7942eSJames Wright   PetscFunctionBeginUser;
24362b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
24462b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
24562b7942eSJames Wright       B[j * N + i] = A[i * M + j];
24662b7942eSJames Wright     }
24762b7942eSJames Wright   }
248d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
24962b7942eSJames Wright }
25062b7942eSJames Wright 
25162b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
252cceb3143SJames Wright static PetscErrorCode SgsDDModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDModelContext *psgsdd_ctx) {
25342454adaSJames Wright   SgsDDModelContext sgsdd_ctx;
25462b7942eSJames Wright   PetscInt          num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
25562b7942eSJames Wright   char              file_path[PETSC_MAX_PATH_LEN];
25662b7942eSJames Wright   PetscScalar      *temp;
25762b7942eSJames Wright 
25862b7942eSJames Wright   PetscFunctionBeginUser;
25962b7942eSJames Wright   {
26042454adaSJames Wright     SgsDDModelContext sgsdd_temp;
26162b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
26262b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
26362b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
26462b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
26562b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
26662b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
26762b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
26862b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
26962b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
27062b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
27162b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
27262b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
27362b7942eSJames Wright   }
27462b7942eSJames Wright 
27562b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
27642454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
27762b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
27842454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
27962b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
28042454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
28162b7942eSJames Wright 
28262b7942eSJames Wright   {
28362b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
28462b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
28542454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
28662b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
28762b7942eSJames Wright     PetscCall(PetscFree(temp));
28862b7942eSJames Wright   }
28962b7942eSJames Wright   {
29062b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
29162b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
29242454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
29362b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
29462b7942eSJames Wright     PetscCall(PetscFree(temp));
29562b7942eSJames Wright   }
29662b7942eSJames Wright 
29762b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
29862b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
299d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
30062b7942eSJames Wright }
30162b7942eSJames Wright 
30242454adaSJames Wright PetscErrorCode SgsDDModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
303ee1455b7SJames Wright   PetscReal                alpha = 0;
30442454adaSJames Wright   SgsDDModelContext        sgsdd_ctx;
30562b7942eSJames Wright   MPI_Comm                 comm                           = user->comm;
306ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
30742454adaSJames Wright   SgsDDModelSetupData      sgs_dd_setup_data;
308ee1455b7SJames Wright   NewtonianIdealGasContext gas;
30962b7942eSJames Wright 
31006f41313SJames Wright   PetscFunctionBeginUser;
311f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, user->phys->state_var, ceed_data->elem_restr_q, ceed_data->basis_q,
312f6ac214eSJames Wright                                             &user->grad_velo_proj));
3139ab09d51SJames Wright 
31462b7942eSJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
31562b7942eSJames Wright 
3164b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
31762b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
31862b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
31962b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
32062b7942eSJames Wright   PetscOptionsEnd();
32162b7942eSJames Wright 
3224b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
32362b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
32462b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
32562b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
32662b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
32762b7942eSJames Wright 
32842454adaSJames Wright   PetscCall(SgsDDModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
32962b7942eSJames Wright 
330ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
331ee1455b7SJames Wright   PetscCall(PetscNew(&user->sgs_dd_data));
33267263decSKenneth E. Jansen   PetscCall(
33342454adaSJames Wright       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));
334ee1455b7SJames Wright 
335c38c977aSJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
336c38c977aSJames Wright 
337b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
338ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
339b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
340b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
341b4c37c5cSJames Wright   PetscCallCeed(ceed,
342b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
343b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
344ee1455b7SJames Wright 
345*40816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfunction_context, &sgs_dd_setup_data->ifunction_qfctx));
346*40816385SJames Wright 
347c38c977aSJames Wright   // -- Compute and store anisotropy tensor
348c38c977aSJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
349c38c977aSJames Wright                                                      &sgs_dd_setup_data->grid_aniso_ceed));
350c38c977aSJames Wright 
351ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
352cceb3143SJames Wright   PetscCall(SgsDDModelSetupNodalEvaluation_Fused(ceed, user, ceed_data, sgs_dd_setup_data));
353ee1455b7SJames Wright 
3549c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
35542454adaSJames Wright   PetscCall(SgsModelSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data));
3569c678832SJames Wright 
35742454adaSJames Wright   PetscCall(SgsDDModelSetupDataDestroy(sgs_dd_setup_data));
358d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
35962b7942eSJames Wright }
360ee1455b7SJames Wright 
36142454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
362ee1455b7SJames Wright   PetscFunctionBeginUser;
363d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
364b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
365ee1455b7SJames Wright 
366b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
367ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
36841edf198SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
369ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
370ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
371d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
372ee1455b7SJames Wright }
373