xref: /honee/problems/sgs_dd_model.c (revision e07531f7cad484157c281bb3cb12a29a67a96955)
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 
11c38c977aSJames Wright typedef struct {
12c38c977aSJames Wright   CeedElemRestriction      elem_restr_grid_aniso, elem_restr_sgs;
13c38c977aSJames Wright   CeedVector               grid_aniso_ceed;
1440816385SJames Wright   CeedQFunctionContext     sgsdd_qfctx, ifunction_qfctx;
154c07ec22SJames Wright   SGSModelDDImplementation sgs_dd_model_implementation;
16ad494f68SJames Wright } *SgsDDSetupData;
17c38c977aSJames Wright 
18ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
19b4c37c5cSJames Wright   Ceed ceed;
20ad494f68SJames Wright 
21c38c977aSJames Wright   PetscFunctionBeginUser;
22b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
23ad494f68SJames Wright 
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));
2840816385SJames 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
34ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(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 
56ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method
57ad494f68SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
58cceb3143SJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
59cceb3143SJames Wright   PetscMemType q_mem_type;
60cceb3143SJames Wright 
61cceb3143SJames Wright   PetscFunctionBeginUser;
62a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
63cceb3143SJames Wright 
64cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
65cceb3143SJames Wright 
66a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
67cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
68cceb3143SJames Wright }
69cceb3143SJames Wright 
70b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
71ad494f68SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
7242454adaSJames Wright   SgsDDData           sgs_dd_data = user->sgs_dd_data;
735930f037SJames Wright   CeedQFunction       qf_sgs_dd_nodal;
745930f037SJames Wright   CeedOperator        op_sgs_dd_nodal;
7515c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
76defe8520SJames Wright   PetscInt            dim;
775930f037SJames Wright   CeedVector          inv_multiplicity;
78ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
7915c18037SJames Wright   DMLabel             domain_label = NULL;
8015c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
81ee1455b7SJames Wright 
82ee1455b7SJames Wright   PetscFunctionBeginUser;
83ee1455b7SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
84b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
85b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
86b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
87ee1455b7SJames Wright 
88ee1455b7SJames Wright   {  // Get velocity gradient information
89ee1455b7SJames Wright     CeedOperatorField op_field;
90b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
91b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
92b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
93ee1455b7SJames Wright   }
9415c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
95b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
96ee1455b7SJames Wright 
975930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
985930f037SJames Wright                                    &inv_multiplicity));
99ee1455b7SJames Wright 
100ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
101ee1455b7SJames Wright   switch (user->phys->state_var) {
102ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
103ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
104ee1455b7SJames Wright       break;
105ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
106ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
107ee1455b7SJames Wright       break;
1089b103f75SJames Wright     case STATEVAR_ENTROPY:
1099b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
1109b103f75SJames Wright       break;
111ee1455b7SJames Wright   }
112ee1455b7SJames Wright 
113ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
114ee1455b7SJames Wright   CeedBasis basis_x_to_q;
115b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
116ee1455b7SJames Wright 
117b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
118b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
119b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
120b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
121b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
122b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
123b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
124ee1455b7SJames Wright 
125b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
12658e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
127b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord));
12858e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
12958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
130b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
13158e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
13258e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
133ee1455b7SJames Wright 
1344b0f6111SJames 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,
1354b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
136ee1455b7SJames Wright 
137ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
138ad494f68SJames Wright   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
139ee1455b7SJames Wright 
140b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
141b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
142b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
143b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
144b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
145d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
146ee1455b7SJames Wright }
147ee1455b7SJames Wright 
1484c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation
1494c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1504c07ec22SJames Wright                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
151b87d60b3SJames Wright                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
152b87d60b3SJames Wright                                                                 void **ctx) {
153b87d60b3SJames Wright   CeedQFunction         qf_sgs_dd_inference;
154b87d60b3SJames Wright   CeedOperator          op_sgs_dd_inference;
155b87d60b3SJames Wright   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
156b87d60b3SJames Wright 
157b87d60b3SJames Wright   PetscFunctionBeginUser;
158b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
159b87d60b3SJames Wright                                                   &qf_sgs_dd_inference));
160b87d60b3SJames Wright 
161b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
162b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
163b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
164b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
165b87d60b3SJames Wright 
166b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
167b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
168b87d60b3SJames Wright   PetscCallCeed(ceed,
169b87d60b3SJames Wright                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
170b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
171b87d60b3SJames Wright 
172b87d60b3SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
173b87d60b3SJames Wright                                        op_context));
174b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
175b87d60b3SJames Wright 
176b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
177b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
178b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
179b87d60b3SJames Wright }
180b87d60b3SJames Wright 
1814c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation
1824c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
183b87d60b3SJames Wright   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
184b87d60b3SJames Wright 
185b87d60b3SJames Wright   PetscFunctionBeginUser;
186b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
187b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
188b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeBegin());
189b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
190b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeEnd());
191b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
192b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
193b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
194b87d60b3SJames Wright }
195b87d60b3SJames Wright 
1964c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch
1974c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1984c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
1994c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
2004c07ec22SJames Wright                                                                  void **ctx) {
2014c07ec22SJames Wright   const char     *ceed_resource;
2024c07ec22SJames Wright   char            model_path[PETSC_MAX_PATH_LEN] = "";
2034c07ec22SJames Wright   TorchDeviceType model_device_type;
2044c07ec22SJames Wright 
2054c07ec22SJames Wright   PetscFunctionBeginUser;
2064c07ec22SJames Wright   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
2074c07ec22SJames Wright   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
2084c07ec22SJames Wright   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
2097ffa0ff8SJames Wright   // On-device XPU is not working reliably currently, default to CPU inference evaluation
2107ffa0ff8SJames Wright   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
2114c07ec22SJames Wright   else model_device_type = TORCH_DEVICE_CPU;
2124c07ec22SJames Wright   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
2134c07ec22SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
2144c07ec22SJames Wright 
2154c07ec22SJames Wright   PetscCall(LoadModel_Torch(model_path, model_device_type));
2164c07ec22SJames Wright 
2174c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2184c07ec22SJames Wright }
2194c07ec22SJames Wright 
2204c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch
2214c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
2224c07ec22SJames Wright   static PetscBool run_through = PETSC_FALSE;
2234c07ec22SJames Wright   PetscFunctionBeginUser;
2244c07ec22SJames Wright   if (!run_through) {
2254c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
2264c07ec22SJames Wright   }
2274c07ec22SJames Wright   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
2284c07ec22SJames Wright   if (!run_through) {
2294c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
2304c07ec22SJames Wright     run_through = PETSC_TRUE;
2314c07ec22SJames Wright   }
2324c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2334c07ec22SJames Wright }
2344c07ec22SJames Wright 
235b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method
236b87d60b3SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
237b87d60b3SJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
238b87d60b3SJames Wright   PetscMemType q_mem_type;
239b87d60b3SJames Wright   Vec          DD_Inputs_loc, DD_Outputs_loc;
240b87d60b3SJames Wright 
241b87d60b3SJames Wright   PetscFunctionBeginUser;
242b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
243b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
244a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
245b87d60b3SJames Wright 
246b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
247b87d60b3SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
248b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
249b87d60b3SJames Wright 
250a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
251b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
252b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
253b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
254b87d60b3SJames Wright }
255b87d60b3SJames Wright 
256b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
257b87d60b3SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
258b87d60b3SJames Wright   SgsDDData           sgs_dd_data = user->sgs_dd_data;
259b87d60b3SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
260b87d60b3SJames Wright   PetscInt            dim;
261b87d60b3SJames Wright   CeedVector          inv_multiplicity, eigvec;
262b87d60b3SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
263b87d60b3SJames Wright       elem_restr_dd_outputs;
264b87d60b3SJames Wright   DMLabel  domain_label = NULL;
265b87d60b3SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
266b87d60b3SJames Wright 
267b87d60b3SJames Wright   PetscFunctionBeginUser;
268b87d60b3SJames Wright   {  // Create DMs for data-driven input and output values
269b87d60b3SJames Wright     PetscSection section;
270b87d60b3SJames Wright     PetscInt     degree, q_extra;
271b87d60b3SJames Wright     {  // Get degree and number of quadrature points from dm_sgs
272b87d60b3SJames Wright       PetscFE         fe;
273b87d60b3SJames Wright       PetscSpace      basis;
274b87d60b3SJames Wright       PetscQuadrature quadrature;
275b87d60b3SJames Wright       PetscInt        num_qpnts;
276b87d60b3SJames Wright       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
277b87d60b3SJames Wright       PetscCall(PetscFEGetBasisSpace(fe, &basis));
278b87d60b3SJames Wright       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
279b87d60b3SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
280b87d60b3SJames Wright       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
281b87d60b3SJames Wright       q_extra = degree - num_qpnts;
282b87d60b3SJames Wright     }
283b87d60b3SJames Wright 
284b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
285b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
286b87d60b3SJames 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));
287b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
288b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
289b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
290b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
291b87d60b3SJames Wright 
292b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
293b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
294b87d60b3SJames Wright     }
295b87d60b3SJames Wright 
296b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
297b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
298b87d60b3SJames 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));
299b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
300b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
301b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
302b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
303b87d60b3SJames Wright 
304b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
305b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
306b87d60b3SJames Wright     }
307b87d60b3SJames Wright   }
308b87d60b3SJames Wright 
309b87d60b3SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
310b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
311b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
312b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
313b87d60b3SJames Wright 
314b87d60b3SJames Wright   {  // Get velocity gradient information
315b87d60b3SJames Wright     CeedOperatorField op_field;
316b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
317b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
318b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
319b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
320b87d60b3SJames Wright   }
321b87d60b3SJames Wright 
322b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
323b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
324b87d60b3SJames Wright   PetscCall(
325b87d60b3SJames Wright       DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec));
326b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
327b87d60b3SJames Wright 
328b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs));
329b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs));
330b87d60b3SJames Wright 
3315930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
3325930f037SJames Wright                                    &inv_multiplicity));
333b87d60b3SJames Wright 
334b87d60b3SJames Wright   {  // Create operator for data-driven input evaluation
335b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_inputs;
336b87d60b3SJames Wright     CeedOperator  op_sgs_dd_inputs;
337b87d60b3SJames Wright 
338b87d60b3SJames Wright     switch (user->phys->state_var) {
339b87d60b3SJames Wright       case STATEVAR_PRIMITIVE:
340b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
341b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
342b87d60b3SJames Wright         break;
343b87d60b3SJames Wright       case STATEVAR_CONSERVATIVE:
344b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
345b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
346b87d60b3SJames Wright         break;
3479b103f75SJames Wright       case STATEVAR_ENTROPY:
3489b103f75SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
3499b103f75SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
3509b103f75SJames Wright         break;
351b87d60b3SJames Wright     }
352b87d60b3SJames Wright 
353b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
354b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
355b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
356b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
357b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
358b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
359b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
360b87d60b3SJames Wright 
361b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
362b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
363b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
364b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
365b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
366b87d60b3SJames Wright     PetscCallCeed(ceed,
367b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
368b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
369b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
370b87d60b3SJames Wright 
371b87d60b3SJames Wright     PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
372b87d60b3SJames Wright                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
373b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
374b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
375b87d60b3SJames Wright   }
376b87d60b3SJames Wright 
377b87d60b3SJames Wright   {  // Create operator for data-driven output handling
378b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_outputs;
379b87d60b3SJames Wright     CeedOperator  op_sgs_dd_outputs;
380b87d60b3SJames Wright 
381b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
382b87d60b3SJames Wright                                                     &qf_sgs_dd_outputs));
383b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
384b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
385b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
386b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
387b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
388b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
389b87d60b3SJames Wright 
390b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
391b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
392b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
393b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
394b87d60b3SJames Wright     PetscCallCeed(ceed,
395b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
396b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
397b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
398b87d60b3SJames Wright 
399b87d60b3SJames 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,
400b87d60b3SJames Wright                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
401b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
402b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
403b87d60b3SJames Wright   }
404b87d60b3SJames Wright 
405b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
4064c07ec22SJames Wright 
4074c07ec22SJames Wright   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
4084c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
4094c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
410b87d60b3SJames Wright                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4114c07ec22SJames Wright   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
4124c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
4134c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
4144c07ec22SJames Wright                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4154c07ec22SJames Wright   }
416b87d60b3SJames Wright 
417b87d60b3SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
418b87d60b3SJames Wright 
419b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
420b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
421b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
422b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
423b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
424b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
425b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
426b87d60b3SJames Wright }
427b87d60b3SJames Wright 
4289c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
429ad494f68SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
43042454adaSJames Wright   SgsDDData     sgs_dd_data = user->sgs_dd_data;
43167263decSKenneth E. Jansen   CeedInt       num_comp_q, num_comp_qd, num_comp_x;
432defe8520SJames Wright   PetscInt      dim;
4339c678832SJames Wright   CeedQFunction qf_sgs_apply;
4349c678832SJames Wright   CeedOperator  op_sgs_apply;
4359c678832SJames Wright   CeedBasis     basis_sgs;
4369c678832SJames Wright 
4379c678832SJames Wright   PetscFunctionBeginUser;
4389c678832SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
439b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
440b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd));
441b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
4429c678832SJames Wright 
44367263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs));
4449c678832SJames Wright 
4459c678832SJames Wright   switch (user->phys->state_var) {
4469c678832SJames Wright     case STATEVAR_PRIMITIVE:
44742454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
4489c678832SJames Wright       break;
4499c678832SJames Wright     case STATEVAR_CONSERVATIVE:
45042454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
4519c678832SJames Wright       break;
4529b103f75SJames Wright     case STATEVAR_ENTROPY:
4539b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
4549b103f75SJames Wright       break;
4559c678832SJames Wright   }
4569c678832SJames Wright 
45740816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
458b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
459b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE));
460b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
461b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
4629c678832SJames Wright 
463b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
464b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
46558e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
466b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
467b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
4689c678832SJames Wright 
4699c678832SJames Wright   PetscCall(
4709c678832SJames 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));
4719c678832SJames Wright 
47287edb941SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs));
473b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
474b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
475d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4769c678832SJames Wright }
4779c678832SJames Wright 
4789c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
479ad494f68SJames Wright PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc) {
48042454adaSJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
4819c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
482cceb3143SJames Wright   PetscMemType sgs_nodal_mem_type;
4839c678832SJames Wright 
4849c678832SJames Wright   PetscFunctionBeginUser;
485b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
4869c678832SJames Wright   PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
487f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionApply(user->grad_velo_proj, Q_loc, VelocityGradient));
4889c678832SJames Wright 
4899c678832SJames Wright   // -- Compute Nodal SGS tensor
4909c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
491cceb3143SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(user, Q_loc, VelocityGradient, SGSNodal_loc));
4929c678832SJames Wright 
4939c678832SJames Wright   // -- Compute contribution of the SGS stress
494a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
4959c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
4969c678832SJames Wright 
4979c678832SJames Wright   // -- Return local SGS vector
498a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
4999c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5009c678832SJames Wright   PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
501b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
502d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5039c678832SJames Wright }
5049c678832SJames Wright 
50562b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
506cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
50762b7942eSJames Wright   PetscFunctionBeginUser;
50862b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
50962b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
51062b7942eSJames Wright       B[j * N + i] = A[i * M + j];
51162b7942eSJames Wright     }
51262b7942eSJames Wright   }
513d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
51462b7942eSJames Wright }
51562b7942eSJames Wright 
51662b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
517ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
518ad494f68SJames Wright   SgsDDContext sgsdd_ctx;
51962b7942eSJames Wright   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
52062b7942eSJames Wright   char         file_path[PETSC_MAX_PATH_LEN];
52162b7942eSJames Wright   PetscScalar *temp;
52262b7942eSJames Wright 
52362b7942eSJames Wright   PetscFunctionBeginUser;
52462b7942eSJames Wright   {
525ad494f68SJames Wright     SgsDDContext sgsdd_temp;
52662b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
52762b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
52862b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
52962b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
53062b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
53162b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
53262b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
53362b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
53462b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
53562b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
53662b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
53762b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
53862b7942eSJames Wright   }
53962b7942eSJames Wright 
54062b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
54142454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
54262b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
54342454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
54462b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
54542454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
54662b7942eSJames Wright 
54762b7942eSJames Wright   {
54862b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
54962b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
55042454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
55162b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
55262b7942eSJames Wright     PetscCall(PetscFree(temp));
55362b7942eSJames Wright   }
55462b7942eSJames Wright   {
55562b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
55662b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
55742454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
55862b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
55962b7942eSJames Wright     PetscCall(PetscFree(temp));
56062b7942eSJames Wright   }
56162b7942eSJames Wright 
56262b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
56362b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
564d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
56562b7942eSJames Wright }
56662b7942eSJames Wright 
567991aef52SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
568ee1455b7SJames Wright   PetscReal                alpha = 0;
569ad494f68SJames Wright   SgsDDContext             sgsdd_ctx;
57062b7942eSJames Wright   MPI_Comm                 comm                           = user->comm;
571ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
572ad494f68SJames Wright   SgsDDSetupData           sgs_dd_setup_data;
573ee1455b7SJames Wright   NewtonianIdealGasContext gas;
57462b7942eSJames Wright 
57506f41313SJames Wright   PetscFunctionBeginUser;
576f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, user->phys->state_var, ceed_data->elem_restr_q, ceed_data->basis_q,
577f6ac214eSJames Wright                                             &user->grad_velo_proj));
5789ab09d51SJames Wright 
579b87d60b3SJames Wright   PetscCall(PetscNew(&user->sgs_dd_data));
580b87d60b3SJames Wright   user->sgs_dd_data->num_comp_inputs  = 6;
581b87d60b3SJames Wright   user->sgs_dd_data->num_comp_outputs = 6;
58262b7942eSJames Wright 
5834c07ec22SJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
5844c07ec22SJames Wright 
5854b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
58662b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
58762b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
58862b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
5894c07ec22SJames Wright   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
5904c07ec22SJames Wright   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
5914c07ec22SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
5924c07ec22SJames Wright                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
5934c07ec22SJames Wright                              NULL));
59462b7942eSJames Wright   PetscOptionsEnd();
59562b7942eSJames Wright 
596b87d60b3SJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
5974b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
59862b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
59962b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
60062b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
60162b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
60262b7942eSJames Wright 
603ad494f68SJames Wright   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
60462b7942eSJames Wright 
605ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
606ad494f68SJames Wright   PetscCall(SgsDDCreateDM(user->dm, &user->sgs_dd_data->dm_sgs, user->app_ctx->degree, user->app_ctx->q_extra, &user->sgs_dd_data->num_comp_sgs));
607ee1455b7SJames Wright 
608*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
609ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
610*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
611b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
612b4c37c5cSJames Wright   PetscCallCeed(ceed,
613b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
614b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
615ee1455b7SJames Wright 
616*e07531f7SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx));
61740816385SJames Wright 
618c38c977aSJames Wright   // -- Compute and store anisotropy tensor
619c38c977aSJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
620c38c977aSJames Wright                                                      &sgs_dd_setup_data->grid_aniso_ceed));
621c38c977aSJames Wright 
622ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
6234c07ec22SJames Wright   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
6244c07ec22SJames Wright     case SGS_MODEL_DD_FUSED:
6254c07ec22SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, user, ceed_data, sgs_dd_setup_data));
6264c07ec22SJames Wright       break;
6274c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_CEED:
6284c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_TORCH:
6294c07ec22SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, user, ceed_data, sgs_dd_setup_data));
6304c07ec22SJames Wright       break;
6314c07ec22SJames Wright   }
632ee1455b7SJames Wright 
6339c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
634ad494f68SJames Wright   PetscCall(SgsSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data));
6359c678832SJames Wright 
636ad494f68SJames Wright   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
637d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
63862b7942eSJames Wright }
639ee1455b7SJames Wright 
64042454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
641ee1455b7SJames Wright   PetscFunctionBeginUser;
642d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
643b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
644ee1455b7SJames Wright 
645b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
646b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed));
647ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
64841edf198SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
649b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx));
650b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx));
651ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
652b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs));
653b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs));
654b87d60b3SJames 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));
655ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
656d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
657ee1455b7SJames Wright }
658