xref: /honee/problems/sgs_dd_model.c (revision cf8f12c89a11a3279f97d0071ccf240f989bc25c)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
362b7942eSJames Wright 
462b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h"
562b7942eSJames Wright 
662b7942eSJames Wright #include <petscdmplex.h>
762b7942eSJames Wright 
8149fb536SJames Wright #include <navierstokes.h>
94c07ec22SJames Wright #include <sgs_model_torch.h>
1062b7942eSJames Wright 
1182baf964SJames Wright typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
1282baf964SJames Wright typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
1382baf964SJames Wright typedef struct {
1482baf964SJames Wright   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
1582baf964SJames Wright   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
1682baf964SJames Wright   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
1782baf964SJames Wright   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
1882baf964SJames Wright   SgsDDNodalStressEval      sgs_nodal_eval;
1982baf964SJames Wright   SgsDDNodalStressInference sgs_nodal_inference;
2082baf964SJames Wright   void                     *sgs_nodal_inference_ctx;
2182baf964SJames Wright   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
2282baf964SJames Wright } *SgsDDData;
2382baf964SJames Wright 
2482baf964SJames Wright // @brief Destroy `SgsDDData` object
2582baf964SJames Wright static PetscErrorCode SgsDDDataDestroy(SgsDDData *sgs_dd_data) {
2682baf964SJames Wright   SgsDDData sgs_dd_data_ = *sgs_dd_data;
2714bd2a07SJames Wright 
2882baf964SJames Wright   PetscFunctionBeginUser;
2982baf964SJames Wright   if (!sgs_dd_data_) PetscFunctionReturn(PETSC_SUCCESS);
3082baf964SJames Wright   Ceed ceed = sgs_dd_data_->op_sgs_apply_ctx->ceed;
3182baf964SJames Wright 
3282baf964SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->sgs_nodal_ceed));
3382baf964SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->grad_velo_ceed));
3482baf964SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_evaluation_ctx));
3582baf964SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_sgs_apply_ctx));
3682baf964SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_inputs_ctx));
3782baf964SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_outputs_ctx));
3882baf964SJames Wright   PetscCall(DMDestroy(&sgs_dd_data_->dm_sgs));
3982baf964SJames Wright   PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_inputs));
4082baf964SJames Wright   PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_outputs));
4182baf964SJames Wright   if (sgs_dd_data_->sgs_nodal_inference_ctx) PetscCall(sgs_dd_data_->sgs_nodal_inference_ctx_destroy(sgs_dd_data_->sgs_nodal_inference_ctx));
4282baf964SJames Wright   PetscCall(PetscFree(sgs_dd_data_));
4382baf964SJames Wright   *sgs_dd_data = NULL;
4482baf964SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4582baf964SJames Wright }
4682baf964SJames Wright 
47c38c977aSJames Wright typedef struct {
48c38c977aSJames Wright   CeedElemRestriction      elem_restr_grid_aniso, elem_restr_sgs;
49c38c977aSJames Wright   CeedVector               grid_aniso_ceed;
5040816385SJames Wright   CeedQFunctionContext     sgsdd_qfctx, ifunction_qfctx;
514c07ec22SJames Wright   SGSModelDDImplementation sgs_dd_model_implementation;
52ad494f68SJames Wright } *SgsDDSetupData;
53c38c977aSJames Wright 
548340219bSJames Wright #define GRAD_VELO_PROJ_KEY "Gradient of Velocity Projection"
5582baf964SJames Wright #define SGS_DD_DATA_KEY "SGS Data Driven Data"
568340219bSJames Wright 
57ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
58b4c37c5cSJames Wright   Ceed ceed;
59ad494f68SJames Wright 
60c38c977aSJames Wright   PetscFunctionBeginUser;
61b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
62ad494f68SJames Wright 
63b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
64b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
65b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
66b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
6740816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
68c38c977aSJames Wright   PetscCall(PetscFree(sgs_dd_setup_data));
69519781aeSJames Wright   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
70d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
71c38c977aSJames Wright }
72c38c977aSJames Wright 
73ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
74ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
75ee1455b7SJames Wright   PetscSection section;
76ee1455b7SJames Wright 
77ee1455b7SJames Wright   PetscFunctionBeginUser;
78ee1455b7SJames Wright   *num_components = 6;
79ee1455b7SJames Wright 
80ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
810dee9b8eSJames Wright   PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE));
82ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
83ee1455b7SJames Wright 
84da4ca0cfSJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
85ee1455b7SJames Wright 
86ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
87ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
88ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
89ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
90ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
91ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
92ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
93ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
94d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
95ee1455b7SJames Wright };
96ee1455b7SJames Wright 
97ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method
980c373b74SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
9982baf964SJames Wright   SgsDDData    sgs_dd_data;
100cceb3143SJames Wright   PetscMemType q_mem_type;
101cceb3143SJames Wright 
102cceb3143SJames Wright   PetscFunctionBeginUser;
1030c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
1040c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
105cceb3143SJames Wright 
106cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
107cceb3143SJames Wright 
1080c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
109cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
110cceb3143SJames Wright }
111cceb3143SJames Wright 
112b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
113e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
11482baf964SJames Wright   SgsDDData           sgs_dd_data;
1155930f037SJames Wright   CeedQFunction       qf_sgs_dd_nodal;
1165930f037SJames Wright   CeedOperator        op_sgs_dd_nodal;
11715c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
118defe8520SJames Wright   PetscInt            dim;
1195930f037SJames Wright   CeedVector          inv_multiplicity;
120*cf8f12c8SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_q;
121e3db12f8SJames Wright   PetscInt            height = 0, dm_field = 0;
1228340219bSJames Wright   NodalProjectionData grad_velo_proj;
123ee1455b7SJames Wright 
124ee1455b7SJames Wright   PetscFunctionBeginUser;
1250c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
126*cf8f12c8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
127e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
128*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
129b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
1300c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
1310c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj));
132ee1455b7SJames Wright 
133ee1455b7SJames Wright   {  // Get velocity gradient information
134ee1455b7SJames Wright     CeedOperatorField op_field;
1358340219bSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
136b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
137b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
138ee1455b7SJames Wright   }
139e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_sgs));
140b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
141ee1455b7SJames Wright 
142e3db12f8SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, PETSC_FALSE,
143e3db12f8SJames Wright                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
144ee1455b7SJames Wright 
145ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
1460c373b74SJames Wright   switch (honee->phys->state_var) {
147ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
148ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
149ee1455b7SJames Wright       break;
150ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
151ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
152ee1455b7SJames Wright       break;
1539b103f75SJames Wright     case STATEVAR_ENTROPY:
1549b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
1559b103f75SJames Wright       break;
156ee1455b7SJames Wright   }
157ee1455b7SJames Wright 
158ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
159*cf8f12c8SJames Wright   CeedBasis basis_x_to_q, basis_q;
160*cf8f12c8SJames Wright   PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
161*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, basis_q, &basis_x_to_q));
162ee1455b7SJames Wright 
163b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
164b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
165b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
166b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
167b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
168b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
169b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
170ee1455b7SJames Wright 
171b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
172*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
173e3663b90SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord));
17458e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
17558e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
176b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
17758e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
17858e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
179ee1455b7SJames Wright 
1808340219bSJames Wright   PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, NULL, sgs_dd_data->sgs_nodal_ceed, NULL, NULL,
1818340219bSJames Wright                                        &sgs_dd_data->op_nodal_evaluation_ctx));
182ee1455b7SJames Wright 
183ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
184ad494f68SJames Wright   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
185ee1455b7SJames Wright 
186b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
187b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
188*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
189*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
190b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
191fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
192b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
193b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
194d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
195ee1455b7SJames Wright }
196ee1455b7SJames Wright 
1974c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation
1984c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1994c07ec22SJames Wright                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
200b87d60b3SJames Wright                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
201b87d60b3SJames Wright                                                                 void **ctx) {
202b87d60b3SJames Wright   CeedQFunction         qf_sgs_dd_inference;
203b87d60b3SJames Wright   CeedOperator          op_sgs_dd_inference;
204b87d60b3SJames Wright   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
205b87d60b3SJames Wright 
206b87d60b3SJames Wright   PetscFunctionBeginUser;
207b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
208b87d60b3SJames Wright                                                   &qf_sgs_dd_inference));
209b87d60b3SJames Wright 
210b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
211b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
212b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
213b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
214b87d60b3SJames Wright 
215b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
216b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
217b87d60b3SJames Wright   PetscCallCeed(ceed,
218b87d60b3SJames Wright                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
219b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
220b87d60b3SJames Wright 
221b87d60b3SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
222b87d60b3SJames Wright                                        op_context));
223b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
224b87d60b3SJames Wright 
225b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
226b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
227b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
228b87d60b3SJames Wright }
229b87d60b3SJames Wright 
2304c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation
2314c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
232b87d60b3SJames Wright   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
233b87d60b3SJames Wright 
234b87d60b3SJames Wright   PetscFunctionBeginUser;
235ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
236ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
237b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeBegin());
238b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
239b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeEnd());
240ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
241ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
242b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
243b87d60b3SJames Wright }
244b87d60b3SJames Wright 
2454c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch
2464c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
2474c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
2484c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
2494c07ec22SJames Wright                                                                  void **ctx) {
2504c07ec22SJames Wright   const char     *ceed_resource;
2514c07ec22SJames Wright   char            model_path[PETSC_MAX_PATH_LEN] = "";
2524c07ec22SJames Wright   TorchDeviceType model_device_type;
2534c07ec22SJames Wright 
2544c07ec22SJames Wright   PetscFunctionBeginUser;
2554c07ec22SJames Wright   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
2564c07ec22SJames Wright   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
2574c07ec22SJames Wright   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
2587ffa0ff8SJames Wright   // On-device XPU is not working reliably currently, default to CPU inference evaluation
2597ffa0ff8SJames Wright   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
2604c07ec22SJames Wright   else model_device_type = TORCH_DEVICE_CPU;
2614c07ec22SJames Wright   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
2624c07ec22SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
2634c07ec22SJames Wright 
2644c07ec22SJames Wright   PetscCall(LoadModel_Torch(model_path, model_device_type));
2654c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2664c07ec22SJames Wright }
2674c07ec22SJames Wright 
2684c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch
2694c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
2704c07ec22SJames Wright   static PetscBool run_through = PETSC_FALSE;
27114bd2a07SJames Wright 
2724c07ec22SJames Wright   PetscFunctionBeginUser;
2734c07ec22SJames Wright   if (!run_through) {
2744c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
2754c07ec22SJames Wright   }
2764c07ec22SJames Wright   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
2774c07ec22SJames Wright   if (!run_through) {
2784c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
2794c07ec22SJames Wright     run_through = PETSC_TRUE;
2804c07ec22SJames Wright   }
2814c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2824c07ec22SJames Wright }
2834c07ec22SJames Wright 
284b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method
2850c373b74SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
28682baf964SJames Wright   SgsDDData    sgs_dd_data;
287b87d60b3SJames Wright   PetscMemType q_mem_type;
288b87d60b3SJames Wright   Vec          DD_Inputs_loc, DD_Outputs_loc;
289b87d60b3SJames Wright 
290b87d60b3SJames Wright   PetscFunctionBeginUser;
2910c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
292b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
293b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
2940c373b74SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed));  // q_ceed is an implicit input
295b87d60b3SJames Wright 
296b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
297b87d60b3SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
298b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
299b87d60b3SJames Wright 
3000c373b74SJames Wright   PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc));
301b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
302b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
303b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
304b87d60b3SJames Wright }
305b87d60b3SJames Wright 
306b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
307e3663b90SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
30882baf964SJames Wright   SgsDDData           sgs_dd_data;
309b87d60b3SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
310b87d60b3SJames Wright   PetscInt            dim;
311b87d60b3SJames Wright   CeedVector          inv_multiplicity, eigvec;
3128340219bSJames Wright   NodalProjectionData grad_velo_proj;
313b87d60b3SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
314*cf8f12c8SJames Wright       elem_restr_dd_outputs, elem_restr_q;
315e3db12f8SJames Wright   PetscInt height = 0, dm_field = 0;
316b87d60b3SJames Wright 
317b87d60b3SJames Wright   PetscFunctionBeginUser;
3180c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
319b87d60b3SJames Wright   {  // Create DMs for data-driven input and output values
320b87d60b3SJames Wright     PetscSection section;
321b87d60b3SJames Wright     PetscInt     degree, q_extra;
322b87d60b3SJames Wright     {  // Get degree and number of quadrature points from dm_sgs
323b87d60b3SJames Wright       PetscFE         fe;
324b87d60b3SJames Wright       PetscSpace      basis;
325b87d60b3SJames Wright       PetscQuadrature quadrature;
326b87d60b3SJames Wright       PetscInt        num_qpnts;
327b87d60b3SJames Wright       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
328b87d60b3SJames Wright       PetscCall(PetscFEGetBasisSpace(fe, &basis));
329b87d60b3SJames Wright       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
330b87d60b3SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
331b87d60b3SJames Wright       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
332b87d60b3SJames Wright       q_extra = degree - num_qpnts;
333b87d60b3SJames Wright     }
334b87d60b3SJames Wright 
335b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
3360dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE));
337b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
338b87d60b3SJames Wright     PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_inputs, sgs_dd_data->dm_dd_inputs));
339b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
340b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
341b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
342b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
343b87d60b3SJames Wright 
344b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
345b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
346b87d60b3SJames Wright     }
347b87d60b3SJames Wright 
348b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
3490dee9b8eSJames Wright     PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE));
350b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
351b87d60b3SJames Wright     PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_outputs, sgs_dd_data->dm_dd_outputs));
352b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
353b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
354b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
355b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
356b87d60b3SJames Wright 
357b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
358b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
359b87d60b3SJames Wright     }
360b87d60b3SJames Wright   }
361b87d60b3SJames Wright 
3620c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
363*cf8f12c8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
364e3663b90SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x));
365*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
366b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
3670c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj));
368b87d60b3SJames Wright 
369b87d60b3SJames Wright   {  // Get velocity gradient information
370b87d60b3SJames Wright     CeedOperatorField op_field;
3718340219bSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
372b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
373b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
374b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
375b87d60b3SJames Wright   }
376b87d60b3SJames Wright 
377e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &elem_restr_sgs));
378b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
379e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, num_comp_eigvec,
3809eadbee4SJames Wright                                                       &elem_restr_eigvec));
381b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
382b87d60b3SJames Wright 
383e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field,
384e3db12f8SJames Wright                                             &elem_restr_dd_inputs));
385e3db12f8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field,
386e3db12f8SJames Wright                                             &elem_restr_dd_outputs));
387b87d60b3SJames Wright 
388e3db12f8SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, PETSC_FALSE,
389e3db12f8SJames Wright                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
390b87d60b3SJames Wright 
391b87d60b3SJames Wright   {  // Create operator for data-driven input evaluation
392b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_inputs;
393b87d60b3SJames Wright     CeedOperator  op_sgs_dd_inputs;
394b87d60b3SJames Wright 
3950c373b74SJames Wright     switch (honee->phys->state_var) {
396b87d60b3SJames Wright       case STATEVAR_PRIMITIVE:
397b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
398b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
399b87d60b3SJames Wright         break;
400b87d60b3SJames Wright       case STATEVAR_CONSERVATIVE:
401b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
402b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
403b87d60b3SJames Wright         break;
4049b103f75SJames Wright       case STATEVAR_ENTROPY:
4059b103f75SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
4069b103f75SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
4079b103f75SJames Wright         break;
408b87d60b3SJames Wright     }
409b87d60b3SJames Wright 
410b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
411b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
412b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
413b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
414b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
415b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
416b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
417b87d60b3SJames Wright 
418b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
419*cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", elem_restr_q, CEED_BASIS_NONE, honee->q_ceed));
420b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
421b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
422b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
423b87d60b3SJames Wright     PetscCallCeed(ceed,
424b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
425b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
426b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
427b87d60b3SJames Wright 
4288340219bSJames Wright     PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
429b87d60b3SJames Wright                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
430b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
431b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
432b87d60b3SJames Wright   }
433b87d60b3SJames Wright 
434b87d60b3SJames Wright   {  // Create operator for data-driven output handling
435b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_outputs;
436b87d60b3SJames Wright     CeedOperator  op_sgs_dd_outputs;
437b87d60b3SJames Wright 
438b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
439b87d60b3SJames Wright                                                     &qf_sgs_dd_outputs));
440b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
441b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
442b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
443b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
444b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
445b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
446b87d60b3SJames Wright 
447b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
448b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
449b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
450b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
451b87d60b3SJames Wright     PetscCallCeed(ceed,
452b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
453b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
454b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
455b87d60b3SJames Wright 
456b87d60b3SJames Wright     PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_outputs, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_outputs, NULL, sgs_dd_data->sgs_nodal_ceed,
457b87d60b3SJames Wright                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
458b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
459b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
460b87d60b3SJames Wright   }
461b87d60b3SJames Wright 
462b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
4634c07ec22SJames Wright 
4644c07ec22SJames Wright   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
4654c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
4664c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
467b87d60b3SJames Wright                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4684c07ec22SJames Wright   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
4694c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
4704c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
4714c07ec22SJames Wright                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4724c07ec22SJames Wright   }
473b87d60b3SJames Wright 
474b87d60b3SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
475b87d60b3SJames Wright 
476b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
477b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
478b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
479b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
480b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
481b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
482fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
483*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
484b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
485b87d60b3SJames Wright }
486b87d60b3SJames Wright 
4879c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
488e3663b90SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) {
48982baf964SJames Wright   SgsDDData           sgs_dd_data;
490*cf8f12c8SJames Wright   CeedInt             q_data_size;
491*cf8f12c8SJames Wright   PetscInt            dim, num_comp_q, num_comp_x;
4929c678832SJames Wright   CeedQFunction       qf_sgs_apply;
4939c678832SJames Wright   CeedOperator        op_sgs_apply;
494*cf8f12c8SJames Wright   CeedBasis           basis_sgs, basis_q;
495be29160dSJames Wright   CeedVector          q_data;
496*cf8f12c8SJames Wright   CeedElemRestriction elem_restr_qd, elem_restr_q;
4979c678832SJames Wright 
4989c678832SJames Wright   PetscFunctionBeginUser;
4990c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
5000c373b74SJames Wright   PetscCall(DMGetDimension(honee->dm, &dim));
501*cf8f12c8SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
502*cf8f12c8SJames Wright   PetscCall(DMGetFieldNumComps(honee->dm, 0, &num_comp_q));
503*cf8f12c8SJames Wright   PetscCall(DMGetCoordinateNumComps(honee->dm, &num_comp_x));
5049c678832SJames Wright 
505be29160dSJames Wright   {
506e3db12f8SJames Wright     PetscInt height = 0, dm_field = 0;
507be29160dSJames Wright 
508*cf8f12c8SJames Wright     PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
509e3db12f8SJames Wright     PetscCall(DMPlexCeedBasisCreate(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &basis_sgs));
510e3db12f8SJames Wright     PetscCall(QDataGet(ceed, sgs_dd_data->dm_sgs, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord,
511e3db12f8SJames Wright                        &elem_restr_qd, &q_data, &q_data_size));
512be29160dSJames Wright   }
5139c678832SJames Wright 
5140c373b74SJames Wright   switch (honee->phys->state_var) {
5159c678832SJames Wright     case STATEVAR_PRIMITIVE:
51642454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
5179c678832SJames Wright       break;
5189c678832SJames Wright     case STATEVAR_CONSERVATIVE:
51942454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
5209c678832SJames Wright       break;
5219b103f75SJames Wright     case STATEVAR_ENTROPY:
5229b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
5239b103f75SJames Wright       break;
5249c678832SJames Wright   }
5259c678832SJames Wright 
52640816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
527b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
528be29160dSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE));
529b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
530b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
5319c678832SJames Wright 
532b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
533*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
534be29160dSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
535b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
536*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
5379c678832SJames Wright 
5389eadbee4SJames Wright   PetscCall(OperatorApplyContextCreate(honee->dm, honee->dm, ceed, op_sgs_apply, honee->q_ceed, honee->g_ceed, NULL, NULL,
5399eadbee4SJames Wright                                        &sgs_dd_data->op_sgs_apply_ctx));
5409c678832SJames Wright 
541be29160dSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
542be29160dSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
54387edb941SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs));
544b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
545b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
546*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
547*cf8f12c8SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
548d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5499c678832SJames Wright }
5509c678832SJames Wright 
5519c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
5520c373b74SJames Wright PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) {
55382baf964SJames Wright   SgsDDData           sgs_dd_data;
5549c678832SJames Wright   Vec                 VelocityGradient, SGSNodal_loc;
555cceb3143SJames Wright   PetscMemType        sgs_nodal_mem_type;
5568340219bSJames Wright   NodalProjectionData grad_velo_proj;
5579c678832SJames Wright 
5589c678832SJames Wright   PetscFunctionBeginUser;
559ea615d4cSJames Wright   PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL));
5600c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data));
5610c70a8bcSJames Wright   PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj));
5628340219bSJames Wright   PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &VelocityGradient));
5638340219bSJames Wright   PetscCall(VelocityGradientProjectionApply(grad_velo_proj, Q_loc, VelocityGradient));
5649c678832SJames Wright 
5659c678832SJames Wright   // -- Compute Nodal SGS tensor
5669c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5670c373b74SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc));
5689c678832SJames Wright 
5699c678832SJames Wright   // -- Compute contribution of the SGS stress
570a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
5719c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
5729c678832SJames Wright 
5739c678832SJames Wright   // -- Return local SGS vector
574a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
5759c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5768340219bSJames Wright   PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &VelocityGradient));
577ea615d4cSJames Wright   PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL));
578d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5799c678832SJames Wright }
5809c678832SJames Wright 
58162b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
582cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
58362b7942eSJames Wright   PetscFunctionBeginUser;
58462b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
58562b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
58662b7942eSJames Wright       B[j * N + i] = A[i * M + j];
58762b7942eSJames Wright     }
58862b7942eSJames Wright   }
589d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
59062b7942eSJames Wright }
59162b7942eSJames Wright 
59262b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
593ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
594ad494f68SJames Wright   SgsDDContext sgsdd_ctx;
59562b7942eSJames Wright   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
59662b7942eSJames Wright   char         file_path[PETSC_MAX_PATH_LEN];
59762b7942eSJames Wright   PetscScalar *temp;
59862b7942eSJames Wright 
59962b7942eSJames Wright   PetscFunctionBeginUser;
60062b7942eSJames Wright   {
601ad494f68SJames Wright     SgsDDContext sgsdd_temp;
60262b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
60362b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
60462b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
60562b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
60662b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
60762b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
60862b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
60962b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
61062b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
61162b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
61262b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
61362b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
61462b7942eSJames Wright   }
61562b7942eSJames Wright 
61662b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
61742454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
61862b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
61942454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
62062b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
62142454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
62262b7942eSJames Wright 
62362b7942eSJames Wright   {
62462b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
62562b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
62642454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
62762b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
62862b7942eSJames Wright     PetscCall(PetscFree(temp));
62962b7942eSJames Wright   }
63062b7942eSJames Wright   {
63162b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
63262b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
63342454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
63462b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
63562b7942eSJames Wright     PetscCall(PetscFree(temp));
63662b7942eSJames Wright   }
63762b7942eSJames Wright 
63862b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
63962b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
640d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
64162b7942eSJames Wright }
64262b7942eSJames Wright 
643e3663b90SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) {
644ee1455b7SJames Wright   PetscReal                alpha = 0;
645ad494f68SJames Wright   SgsDDContext             sgsdd_ctx;
6460c373b74SJames Wright   MPI_Comm                 comm                           = honee->comm;
647ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
648ad494f68SJames Wright   SgsDDSetupData           sgs_dd_setup_data;
649cde3d787SJames Wright   NewtonianIdealGasContext newt_ctx;
6508340219bSJames Wright   NodalProjectionData      grad_velo_proj;
65182baf964SJames Wright   SgsDDData                sgs_dd_data;
65262b7942eSJames Wright 
65306f41313SJames Wright   PetscFunctionBeginUser;
654*cf8f12c8SJames Wright   {
655*cf8f12c8SJames Wright     CeedElemRestriction elem_restr_q;
656*cf8f12c8SJames Wright     CeedBasis           basis_q;
657*cf8f12c8SJames Wright 
658*cf8f12c8SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
659*cf8f12c8SJames Wright     PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
660*cf8f12c8SJames Wright     // TODO: Should probably move the elem_restr_q and basis_q creation to inside the velocity gradient projection setup???
661*cf8f12c8SJames Wright     PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, elem_restr_q, basis_q, &grad_velo_proj));
662*cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
663*cf8f12c8SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
664*cf8f12c8SJames Wright   }
6650c70a8bcSJames Wright   PetscCall(HoneeSetContainer(honee, GRAD_VELO_PROJ_KEY, grad_velo_proj, (PetscCtxDestroyFn *)NodalProjectionDataDestroy));
6669ab09d51SJames Wright 
66782baf964SJames Wright   PetscCall(PetscNew(&sgs_dd_data));
66882baf964SJames Wright   sgs_dd_data->num_comp_inputs  = 6;
66982baf964SJames Wright   sgs_dd_data->num_comp_outputs = 6;
6700c70a8bcSJames Wright   PetscCall(HoneeSetContainer(honee, SGS_DD_DATA_KEY, sgs_dd_data, (PetscCtxDestroyFn *)SgsDDDataDestroy));
67162b7942eSJames Wright 
6724c07ec22SJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
6734c07ec22SJames Wright 
6744b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
67562b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
67662b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
67762b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
6784c07ec22SJames Wright   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
6794c07ec22SJames Wright   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
6804c07ec22SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
6814c07ec22SJames Wright                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
6824c07ec22SJames Wright                              NULL));
68362b7942eSJames Wright   PetscOptionsEnd();
68462b7942eSJames Wright 
685b87d60b3SJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
686f5dc303cSJames Wright   *sgsdd_ctx = (struct SgsDDContext_){
687f5dc303cSJames Wright       .num_layers  = 1,
688f5dc303cSJames Wright       .num_inputs  = 6,
689f5dc303cSJames Wright       .num_outputs = 6,
690f5dc303cSJames Wright       .num_neurons = 20,
691f5dc303cSJames Wright       .alpha       = alpha,
692f5dc303cSJames Wright   };
69362b7942eSJames Wright 
694ad494f68SJames Wright   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
69562b7942eSJames Wright 
696ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
69782baf964SJames Wright   PetscCall(SgsDDCreateDM(honee->dm, &sgs_dd_data->dm_sgs, honee->app_ctx->degree, honee->app_ctx->q_extra, &sgs_dd_data->num_comp_sgs));
698ee1455b7SJames Wright 
699cde3d787SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx));
700cde3d787SJames Wright   sgsdd_ctx->newt_ctx = *newt_ctx;
701cde3d787SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx));
7020c373b74SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
703b4c37c5cSJames Wright   PetscCallCeed(ceed,
704b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
705b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
706ee1455b7SJames Wright 
707e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx));
70840816385SJames Wright 
709c38c977aSJames Wright   // -- Compute and store anisotropy tensor
710e3663b90SJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed));
711c38c977aSJames Wright 
712ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
7134c07ec22SJames Wright   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
7144c07ec22SJames Wright     case SGS_MODEL_DD_FUSED:
715e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data));
7164c07ec22SJames Wright       break;
7174c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_CEED:
7184c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_TORCH:
719e3663b90SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data));
7204c07ec22SJames Wright       break;
7214c07ec22SJames Wright   }
722ee1455b7SJames Wright 
7239c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
724e3663b90SJames Wright   PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data));
7259c678832SJames Wright 
726ad494f68SJames Wright   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
727d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
72862b7942eSJames Wright }
729