xref: /honee/problems/sgs_dd_model.c (revision 9b103f75867128bb395d4431a2dd4da8eacd1da9)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
262b7942eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
362b7942eSJames Wright //
462b7942eSJames Wright // SPDX-License-Identifier: BSD-2-Clause
562b7942eSJames Wright //
662b7942eSJames Wright // This file is part of CEED:  http://github.com/ceed
762b7942eSJames Wright 
862b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h"
962b7942eSJames Wright 
1062b7942eSJames Wright #include <petscdmplex.h>
1162b7942eSJames Wright 
124c07ec22SJames Wright #include <sgs_model_torch.h>
1362b7942eSJames Wright #include "../navierstokes.h"
1462b7942eSJames Wright 
15c38c977aSJames Wright typedef struct {
16c38c977aSJames Wright   CeedElemRestriction      elem_restr_grid_aniso, elem_restr_sgs;
17c38c977aSJames Wright   CeedVector               grid_aniso_ceed;
1840816385SJames Wright   CeedQFunctionContext     sgsdd_qfctx, ifunction_qfctx;
194c07ec22SJames Wright   SGSModelDDImplementation sgs_dd_model_implementation;
20ad494f68SJames Wright } *SgsDDSetupData;
21c38c977aSJames Wright 
22ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) {
23b4c37c5cSJames Wright   Ceed ceed;
24ad494f68SJames Wright 
25c38c977aSJames Wright   PetscFunctionBeginUser;
26b4c37c5cSJames Wright   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed));
27ad494f68SJames Wright 
28b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso));
29b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs));
30b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed));
31b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx));
3240816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx));
33c38c977aSJames Wright   PetscCall(PetscFree(sgs_dd_setup_data));
34d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
35c38c977aSJames Wright }
36c38c977aSJames Wright 
37ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes
38ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
39ee1455b7SJames Wright   PetscSection section;
40ee1455b7SJames Wright 
41ee1455b7SJames Wright   PetscFunctionBeginUser;
42ee1455b7SJames Wright   *num_components = 6;
43ee1455b7SJames Wright 
44ee1455b7SJames Wright   PetscCall(DMClone(dm_source, dm_sgs));
45ee1455b7SJames Wright   PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection"));
46ee1455b7SJames Wright 
47da4ca0cfSJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs));
48ee1455b7SJames Wright 
49ee1455b7SJames Wright   PetscCall(DMGetLocalSection(*dm_sgs, &section));
50ee1455b7SJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
51ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX"));
52ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY"));
53ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ"));
54ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ"));
55ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ"));
56ee1455b7SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY"));
57d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
58ee1455b7SJames Wright };
59ee1455b7SJames Wright 
60ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method
61ad494f68SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
62cceb3143SJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
63cceb3143SJames Wright   PetscMemType q_mem_type;
64cceb3143SJames Wright 
65cceb3143SJames Wright   PetscFunctionBeginUser;
66a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
67cceb3143SJames Wright 
68cceb3143SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx));
69cceb3143SJames Wright 
70a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
71cceb3143SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
72cceb3143SJames Wright }
73cceb3143SJames Wright 
74b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator
75ad494f68SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
7642454adaSJames Wright   SgsDDData           sgs_dd_data = user->sgs_dd_data;
775930f037SJames Wright   CeedQFunction       qf_sgs_dd_nodal;
785930f037SJames Wright   CeedOperator        op_sgs_dd_nodal;
7915c18037SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso;
80defe8520SJames Wright   PetscInt            dim;
815930f037SJames Wright   CeedVector          inv_multiplicity;
82ee1455b7SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs;
8315c18037SJames Wright   DMLabel             domain_label = NULL;
8415c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
85ee1455b7SJames Wright 
86ee1455b7SJames Wright   PetscFunctionBeginUser;
87ee1455b7SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
88b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
89b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
90b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
91ee1455b7SJames Wright 
92ee1455b7SJames Wright   {  // Get velocity gradient information
93ee1455b7SJames Wright     CeedOperatorField op_field;
94b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
95b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
96b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
97ee1455b7SJames Wright   }
9815c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
99b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
100ee1455b7SJames Wright 
1015930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
1025930f037SJames Wright                                    &inv_multiplicity));
103ee1455b7SJames Wright 
104ee1455b7SJames Wright   // -- Create operator for SGS DD model nodal evaluation
105ee1455b7SJames Wright   switch (user->phys->state_var) {
106ee1455b7SJames Wright     case STATEVAR_PRIMITIVE:
107ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal));
108ee1455b7SJames Wright       break;
109ee1455b7SJames Wright     case STATEVAR_CONSERVATIVE:
110ad494f68SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal));
111ee1455b7SJames Wright       break;
112*9b103f75SJames Wright     case STATEVAR_ENTROPY:
113*9b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal));
114*9b103f75SJames Wright       break;
115ee1455b7SJames Wright   }
116ee1455b7SJames Wright 
117ee1455b7SJames Wright   // Mesh/geometry order and solution basis order may differ, therefore must interpolate
118ee1455b7SJames Wright   CeedBasis basis_x_to_q;
119b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q));
120ee1455b7SJames Wright 
121b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx));
122b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE));
123b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP));
124b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
125b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
126b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE));
127b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
128ee1455b7SJames Wright 
129b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal));
13058e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
131b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord));
13258e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
13358e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
134b4c37c5cSJames Wright                                            sgs_dd_setup_data->grid_aniso_ceed));
13558e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
13658e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
137ee1455b7SJames Wright 
1384b0f6111SJames 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,
1394b0f6111SJames Wright                                        NULL, &sgs_dd_data->op_nodal_evaluation_ctx));
140ee1455b7SJames Wright 
141ee1455b7SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
142ad494f68SJames Wright   sgs_dd_data->sgs_nodal_eval       = SgsDDNodalStressEval_Fused;
143ee1455b7SJames Wright 
144b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
145b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q));
146b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
147b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal));
148b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal));
149d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
150ee1455b7SJames Wright }
151ee1455b7SJames Wright 
1524c07ec22SJames Wright // @brief Setup data-driven model inference using libCEED native implementation
1534c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
1544c07ec22SJames Wright                                                                 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
155b87d60b3SJames Wright                                                                 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
156b87d60b3SJames Wright                                                                 void **ctx) {
157b87d60b3SJames Wright   CeedQFunction         qf_sgs_dd_inference;
158b87d60b3SJames Wright   CeedOperator          op_sgs_dd_inference;
159b87d60b3SJames Wright   OperatorApplyContext *op_context = (OperatorApplyContext *)ctx;
160b87d60b3SJames Wright 
161b87d60b3SJames Wright   PetscFunctionBeginUser;
162b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc,
163b87d60b3SJames Wright                                                   &qf_sgs_dd_inference));
164b87d60b3SJames Wright 
165b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx));
166b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
167b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE));
168b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
169b87d60b3SJames Wright 
170b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference));
171b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
172b87d60b3SJames Wright   PetscCallCeed(ceed,
173b87d60b3SJames Wright                 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
174b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
175b87d60b3SJames Wright 
176b87d60b3SJames Wright   PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL,
177b87d60b3SJames Wright                                        op_context));
178b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy;
179b87d60b3SJames Wright 
180b87d60b3SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference));
181b87d60b3SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference));
182b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
183b87d60b3SJames Wright }
184b87d60b3SJames Wright 
1854c07ec22SJames Wright // @brief Perform data-driven model inference using libCEED native implementation
1864c07ec22SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
187b87d60b3SJames Wright   OperatorApplyContext op_context = *(OperatorApplyContext *)ctx;
188b87d60b3SJames Wright 
189b87d60b3SJames Wright   PetscFunctionBeginUser;
190b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
191b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
192b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeBegin());
193b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context));
194b40a7e63SJames Wright   PetscCall(PetscLogGpuTimeEnd());
195b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
196b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL));
197b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
198b87d60b3SJames Wright }
199b87d60b3SJames Wright 
2004c07ec22SJames Wright // @brief Setup data-driven model inference using libtorch
2014c07ec22SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data,
2024c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs,
2034c07ec22SJames Wright                                                                  CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity,
2044c07ec22SJames Wright                                                                  void **ctx) {
2054c07ec22SJames Wright   const char     *ceed_resource;
2064c07ec22SJames Wright   char            model_path[PETSC_MAX_PATH_LEN] = "";
2074c07ec22SJames Wright   TorchDeviceType model_device_type;
2084c07ec22SJames Wright 
2094c07ec22SJames Wright   PetscFunctionBeginUser;
2104c07ec22SJames Wright   PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource));
2114c07ec22SJames Wright   if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA;
2124c07ec22SJames Wright   else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP;
2137ffa0ff8SJames Wright   // On-device XPU is not working reliably currently, default to CPU inference evaluation
2147ffa0ff8SJames Wright   // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU;
2154c07ec22SJames Wright   else model_device_type = TORCH_DEVICE_CPU;
2164c07ec22SJames Wright   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL));
2174c07ec22SJames Wright   PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL));
2184c07ec22SJames Wright 
2194c07ec22SJames Wright   PetscCall(LoadModel_Torch(model_path, model_device_type));
2204c07ec22SJames Wright 
2214c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2224c07ec22SJames Wright }
2234c07ec22SJames Wright 
2244c07ec22SJames Wright // @brief Perform data-driven model inference using libtorch
2254c07ec22SJames Wright static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) {
2264c07ec22SJames Wright   static PetscBool run_through = PETSC_FALSE;
2274c07ec22SJames Wright   PetscFunctionBeginUser;
2284c07ec22SJames Wright   if (!run_through) {
2294c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view"));
2304c07ec22SJames Wright   }
2314c07ec22SJames Wright   PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc));
2324c07ec22SJames Wright   if (!run_through) {
2334c07ec22SJames Wright     PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view"));
2344c07ec22SJames Wright     run_through = PETSC_TRUE;
2354c07ec22SJames Wright   }
2364c07ec22SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2374c07ec22SJames Wright }
2384c07ec22SJames Wright 
239b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method
240b87d60b3SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) {
241b87d60b3SJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
242b87d60b3SJames Wright   PetscMemType q_mem_type;
243b87d60b3SJames Wright   Vec          DD_Inputs_loc, DD_Outputs_loc;
244b87d60b3SJames Wright 
245b87d60b3SJames Wright   PetscFunctionBeginUser;
246b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
247b87d60b3SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
248a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, user->q_ceed));  // q_ceed is an implicit input
249b87d60b3SJames Wright 
250b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx));
251b87d60b3SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx));
252b87d60b3SJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx));
253b87d60b3SJames Wright 
254a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(user->q_ceed, q_mem_type, Q_loc));
255b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc));
256b87d60b3SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc));
257b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
258b87d60b3SJames Wright }
259b87d60b3SJames Wright 
260b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators
261b87d60b3SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
262b87d60b3SJames Wright   SgsDDData           sgs_dd_data = user->sgs_dd_data;
263b87d60b3SJames Wright   CeedInt             num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1;
264b87d60b3SJames Wright   PetscInt            dim;
265b87d60b3SJames Wright   CeedVector          inv_multiplicity, eigvec;
266b87d60b3SJames Wright   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs,
267b87d60b3SJames Wright       elem_restr_dd_outputs;
268b87d60b3SJames Wright   DMLabel  domain_label = NULL;
269b87d60b3SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
270b87d60b3SJames Wright 
271b87d60b3SJames Wright   PetscFunctionBeginUser;
272b87d60b3SJames Wright   {  // Create DMs for data-driven input and output values
273b87d60b3SJames Wright     PetscSection section;
274b87d60b3SJames Wright     PetscInt     degree, q_extra;
275b87d60b3SJames Wright     {  // Get degree and number of quadrature points from dm_sgs
276b87d60b3SJames Wright       PetscFE         fe;
277b87d60b3SJames Wright       PetscSpace      basis;
278b87d60b3SJames Wright       PetscQuadrature quadrature;
279b87d60b3SJames Wright       PetscInt        num_qpnts;
280b87d60b3SJames Wright       PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe));
281b87d60b3SJames Wright       PetscCall(PetscFEGetBasisSpace(fe, &basis));
282b87d60b3SJames Wright       PetscCall(PetscSpaceGetDegree(basis, &degree, NULL));
283b87d60b3SJames Wright       PetscCall(PetscFEGetQuadrature(fe, &quadrature));
284b87d60b3SJames Wright       PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts));
285b87d60b3SJames Wright       q_extra = degree - num_qpnts;
286b87d60b3SJames Wright     }
287b87d60b3SJames Wright 
288b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs));
289b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs"));
290b87d60b3SJames 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));
291b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, &section));
292b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
293b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) {
294b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
295b87d60b3SJames Wright 
296b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1));
297b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
298b87d60b3SJames Wright     }
299b87d60b3SJames Wright 
300b87d60b3SJames Wright     PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs));
301b87d60b3SJames Wright     PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs"));
302b87d60b3SJames 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));
303b87d60b3SJames Wright     PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, &section));
304b87d60b3SJames Wright     PetscCall(PetscSectionSetFieldName(section, 0, ""));
305b87d60b3SJames Wright     for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) {
306b87d60b3SJames Wright       char component_name[PETSC_MAX_PATH_LEN];
307b87d60b3SJames Wright 
308b87d60b3SJames Wright       PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1));
309b87d60b3SJames Wright       PetscCall(PetscSectionSetComponentName(section, 0, i, component_name));
310b87d60b3SJames Wright     }
311b87d60b3SJames Wright   }
312b87d60b3SJames Wright 
313b87d60b3SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
314b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
315b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
316b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
317b87d60b3SJames Wright 
318b87d60b3SJames Wright   {  // Get velocity gradient information
319b87d60b3SJames Wright     CeedOperatorField op_field;
320b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
321b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
322b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
323b87d60b3SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL));
324b87d60b3SJames Wright   }
325b87d60b3SJames Wright 
326b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs));
327b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL));
328b87d60b3SJames Wright   PetscCall(
329b87d60b3SJames Wright       DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec));
330b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL));
331b87d60b3SJames Wright 
332b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs));
333b87d60b3SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs));
334b87d60b3SJames Wright 
3355930f037SJames Wright   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity,
3365930f037SJames Wright                                    &inv_multiplicity));
337b87d60b3SJames Wright 
338b87d60b3SJames Wright   {  // Create operator for data-driven input evaluation
339b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_inputs;
340b87d60b3SJames Wright     CeedOperator  op_sgs_dd_inputs;
341b87d60b3SJames Wright 
342b87d60b3SJames Wright     switch (user->phys->state_var) {
343b87d60b3SJames Wright       case STATEVAR_PRIMITIVE:
344b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim,
345b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs));
346b87d60b3SJames Wright         break;
347b87d60b3SJames Wright       case STATEVAR_CONSERVATIVE:
348b87d60b3SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv,
349b87d60b3SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs));
350b87d60b3SJames Wright         break;
351*9b103f75SJames Wright       case STATEVAR_ENTROPY:
352*9b103f75SJames Wright         PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy,
353*9b103f75SJames Wright                                                         ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs));
354*9b103f75SJames Wright         break;
355b87d60b3SJames Wright     }
356b87d60b3SJames Wright 
357b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx));
358b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE));
359b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
360b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
361b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
362b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
363b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE));
364b87d60b3SJames Wright 
365b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs));
366b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed));
367b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
368b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
369b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
370b87d60b3SJames Wright     PetscCallCeed(ceed,
371b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
372b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
373b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
374b87d60b3SJames Wright 
375b87d60b3SJames Wright     PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL,
376b87d60b3SJames Wright                                          &sgs_dd_data->op_nodal_dd_inputs_ctx));
377b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs));
378b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs));
379b87d60b3SJames Wright   }
380b87d60b3SJames Wright 
381b87d60b3SJames Wright   {  // Create operator for data-driven output handling
382b87d60b3SJames Wright     CeedQFunction qf_sgs_dd_outputs;
383b87d60b3SJames Wright     CeedOperator  op_sgs_dd_outputs;
384b87d60b3SJames Wright 
385b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc,
386b87d60b3SJames Wright                                                     &qf_sgs_dd_outputs));
387b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx));
388b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE));
389b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
390b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE));
391b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE));
392b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE));
393b87d60b3SJames Wright 
394b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs));
395b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
396b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
397b87d60b3SJames Wright                                              sgs_dd_setup_data->grid_aniso_ceed));
398b87d60b3SJames Wright     PetscCallCeed(ceed,
399b87d60b3SJames Wright                   CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
400b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec));
401b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
402b87d60b3SJames Wright 
403b87d60b3SJames 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,
404b87d60b3SJames Wright                                          NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx));
405b87d60b3SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs));
406b87d60b3SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs));
407b87d60b3SJames Wright   }
408b87d60b3SJames Wright 
409b87d60b3SJames Wright   sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential;
4104c07ec22SJames Wright 
4114c07ec22SJames Wright   if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) {
4124c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed;
4134c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
414b87d60b3SJames Wright                                                         elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4154c07ec22SJames Wright   } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) {
4164c07ec22SJames Wright     sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch;
4174c07ec22SJames Wright     PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs,
4184c07ec22SJames Wright                                                          elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx));
4194c07ec22SJames Wright   }
420b87d60b3SJames Wright 
421b87d60b3SJames Wright   sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs;
422b87d60b3SJames Wright 
423b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
424b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&eigvec));
425b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
426b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec));
427b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs));
428b87d60b3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs));
429b87d60b3SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
430b87d60b3SJames Wright }
431b87d60b3SJames Wright 
4329c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual
433ad494f68SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) {
43442454adaSJames Wright   SgsDDData     sgs_dd_data = user->sgs_dd_data;
43567263decSKenneth E. Jansen   CeedInt       num_comp_q, num_comp_qd, num_comp_x;
436defe8520SJames Wright   PetscInt      dim;
4379c678832SJames Wright   CeedQFunction qf_sgs_apply;
4389c678832SJames Wright   CeedOperator  op_sgs_apply;
4399c678832SJames Wright   CeedBasis     basis_sgs;
4409c678832SJames Wright 
4419c678832SJames Wright   PetscFunctionBeginUser;
4429c678832SJames Wright   PetscCall(DMGetDimension(user->dm, &dim));
443b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q));
444b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd));
445b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x));
4469c678832SJames Wright 
44767263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs));
4489c678832SJames Wright 
4499c678832SJames Wright   switch (user->phys->state_var) {
4509c678832SJames Wright     case STATEVAR_PRIMITIVE:
45142454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply));
4529c678832SJames Wright       break;
4539c678832SJames Wright     case STATEVAR_CONSERVATIVE:
45442454adaSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply));
4559c678832SJames Wright       break;
456*9b103f75SJames Wright     case STATEVAR_ENTROPY:
457*9b103f75SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply));
458*9b103f75SJames Wright       break;
4599c678832SJames Wright   }
4609c678832SJames Wright 
46140816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx));
462b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP));
463b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE));
464b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP));
465b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
4669c678832SJames Wright 
467b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply));
468b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
46958e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
470b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed));
471b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
4729c678832SJames Wright 
4739c678832SJames Wright   PetscCall(
4749c678832SJames 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));
4759c678832SJames Wright 
476b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply));
477b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply));
478d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
4799c678832SJames Wright }
4809c678832SJames Wright 
4819c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual
482ad494f68SJames Wright PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc) {
48342454adaSJames Wright   SgsDDData    sgs_dd_data = user->sgs_dd_data;
4849c678832SJames Wright   Vec          VelocityGradient, SGSNodal_loc;
485cceb3143SJames Wright   PetscMemType sgs_nodal_mem_type;
4869c678832SJames Wright 
4879c678832SJames Wright   PetscFunctionBeginUser;
488b40a7e63SJames Wright   PetscCall(PetscLogEventBegin(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
4899c678832SJames Wright   PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
490f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionApply(user->grad_velo_proj, Q_loc, VelocityGradient));
4919c678832SJames Wright 
4929c678832SJames Wright   // -- Compute Nodal SGS tensor
4939c678832SJames Wright   PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
494cceb3143SJames Wright   PetscCall(sgs_dd_data->sgs_nodal_eval(user, Q_loc, VelocityGradient, SGSNodal_loc));
4959c678832SJames Wright 
4969c678832SJames Wright   // -- Compute contribution of the SGS stress
497a7dac1d5SJames Wright   PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed));  // sgs_nodal_ceed is an implicit input
4989c678832SJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx));
4999c678832SJames Wright 
5009c678832SJames Wright   // -- Return local SGS vector
501a7dac1d5SJames Wright   PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc));
5029c678832SJames Wright   PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc));
5039c678832SJames Wright   PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient));
504b40a7e63SJames Wright   PetscCall(PetscLogEventEnd(FLUIDS_SgsModel, Q_loc, G_loc, NULL, NULL));
505d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
5069c678832SJames Wright }
5079c678832SJames Wright 
50862b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN
509cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) {
51062b7942eSJames Wright   PetscFunctionBeginUser;
51162b7942eSJames Wright   for (PetscInt i = 0; i < N; i++) {
51262b7942eSJames Wright     for (PetscInt j = 0; j < M; j++) {
51362b7942eSJames Wright       B[j * N + i] = A[i * M + j];
51462b7942eSJames Wright     }
51562b7942eSJames Wright   }
516d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
51762b7942eSJames Wright }
51862b7942eSJames Wright 
51962b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct
520ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) {
521ad494f68SJames Wright   SgsDDContext sgsdd_ctx;
52262b7942eSJames Wright   PetscInt     num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons;
52362b7942eSJames Wright   char         file_path[PETSC_MAX_PATH_LEN];
52462b7942eSJames Wright   PetscScalar *temp;
52562b7942eSJames Wright 
52662b7942eSJames Wright   PetscFunctionBeginUser;
52762b7942eSJames Wright   {
528ad494f68SJames Wright     SgsDDContext sgsdd_temp;
52962b7942eSJames Wright     PetscCall(PetscNew(&sgsdd_temp));
53062b7942eSJames Wright     *sgsdd_temp                     = **psgsdd_ctx;
53162b7942eSJames Wright     sgsdd_temp->offsets.bias1       = 0;
53262b7942eSJames Wright     sgsdd_temp->offsets.bias2       = sgsdd_temp->offsets.bias1 + num_neurons;
53362b7942eSJames Wright     sgsdd_temp->offsets.weight1     = sgsdd_temp->offsets.bias2 + num_neurons;
53462b7942eSJames Wright     sgsdd_temp->offsets.weight2     = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs;
53562b7942eSJames Wright     sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons;
53662b7942eSJames Wright     PetscInt total_num_scalars      = sgsdd_temp->offsets.out_scaling + 2 * num_outputs;
53762b7942eSJames Wright     sgsdd_temp->total_bytes         = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]);
53862b7942eSJames Wright     PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx));
53962b7942eSJames Wright     *sgsdd_ctx = *sgsdd_temp;
54062b7942eSJames Wright     PetscCall(PetscFree(sgsdd_temp));
54162b7942eSJames Wright   }
54262b7942eSJames Wright 
54362b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat"));
54442454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1]));
54562b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat"));
54642454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2]));
54762b7942eSJames Wright   PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat"));
54842454adaSJames Wright   PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling]));
54962b7942eSJames Wright 
55062b7942eSJames Wright   {
55162b7942eSJames Wright     PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp));
55262b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat"));
55342454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
55462b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons));
55562b7942eSJames Wright     PetscCall(PetscFree(temp));
55662b7942eSJames Wright   }
55762b7942eSJames Wright   {
55862b7942eSJames Wright     PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp));
55962b7942eSJames Wright     PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat"));
56042454adaSJames Wright     PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp));
56162b7942eSJames Wright     PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs));
56262b7942eSJames Wright     PetscCall(PetscFree(temp));
56362b7942eSJames Wright   }
56462b7942eSJames Wright 
56562b7942eSJames Wright   PetscCall(PetscFree(*psgsdd_ctx));
56662b7942eSJames Wright   *psgsdd_ctx = sgsdd_ctx;
567d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
56862b7942eSJames Wright }
56962b7942eSJames Wright 
570991aef52SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
571ee1455b7SJames Wright   PetscReal                alpha = 0;
572ad494f68SJames Wright   SgsDDContext             sgsdd_ctx;
57362b7942eSJames Wright   MPI_Comm                 comm                           = user->comm;
574ee1455b7SJames Wright   char                     sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters";
575ad494f68SJames Wright   SgsDDSetupData           sgs_dd_setup_data;
576ee1455b7SJames Wright   NewtonianIdealGasContext gas;
57762b7942eSJames Wright 
57806f41313SJames Wright   PetscFunctionBeginUser;
579f6ac214eSJames Wright   PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, user->phys->state_var, ceed_data->elem_restr_q, ceed_data->basis_q,
580f6ac214eSJames Wright                                             &user->grad_velo_proj));
5819ab09d51SJames Wright 
582b87d60b3SJames Wright   PetscCall(PetscNew(&user->sgs_dd_data));
583b87d60b3SJames Wright   user->sgs_dd_data->num_comp_inputs  = 6;
584b87d60b3SJames Wright   user->sgs_dd_data->num_comp_outputs = 6;
58562b7942eSJames Wright 
5864c07ec22SJames Wright   PetscCall(PetscNew(&sgs_dd_setup_data));
5874c07ec22SJames Wright 
5884b0f6111SJames Wright   PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL);
58962b7942eSJames Wright   PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL));
59062b7942eSJames Wright   PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir,
59162b7942eSJames Wright                                sgs_dd_dir, sizeof(sgs_dd_dir), NULL));
5924c07ec22SJames Wright   PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead"));
5934c07ec22SJames Wright   sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED;
5944c07ec22SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations,
5954c07ec22SJames Wright                              (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation,
5964c07ec22SJames Wright                              NULL));
59762b7942eSJames Wright   PetscOptionsEnd();
59862b7942eSJames Wright 
599b87d60b3SJames Wright   PetscCall(PetscNew(&sgsdd_ctx));
6004b0f6111SJames Wright   sgsdd_ctx->num_layers  = 1;
60162b7942eSJames Wright   sgsdd_ctx->num_inputs  = 6;
60262b7942eSJames Wright   sgsdd_ctx->num_outputs = 6;
60362b7942eSJames Wright   sgsdd_ctx->num_neurons = 20;
60462b7942eSJames Wright   sgsdd_ctx->alpha       = alpha;
60562b7942eSJames Wright 
606ad494f68SJames Wright   PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx));
60762b7942eSJames Wright 
608ee1455b7SJames Wright   // -- Create DM for storing SGS tensor at nodes
609ad494f68SJames 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));
610ee1455b7SJames Wright 
611b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
612ee1455b7SJames Wright   sgsdd_ctx->gas = *gas;
613b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
614b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx));
615b4c37c5cSJames Wright   PetscCallCeed(ceed,
616b4c37c5cSJames Wright                 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx));
617b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc));
618ee1455b7SJames Wright 
61940816385SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfunction_context, &sgs_dd_setup_data->ifunction_qfctx));
62040816385SJames Wright 
621c38c977aSJames Wright   // -- Compute and store anisotropy tensor
622c38c977aSJames Wright   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso,
623c38c977aSJames Wright                                                      &sgs_dd_setup_data->grid_aniso_ceed));
624c38c977aSJames Wright 
625ee1455b7SJames Wright   // -- Create Nodal Evaluation Operator
6264c07ec22SJames Wright   switch (sgs_dd_setup_data->sgs_dd_model_implementation) {
6274c07ec22SJames Wright     case SGS_MODEL_DD_FUSED:
6284c07ec22SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, user, ceed_data, sgs_dd_setup_data));
6294c07ec22SJames Wright       break;
6304c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_CEED:
6314c07ec22SJames Wright     case SGS_MODEL_DD_SEQENTIAL_TORCH:
6324c07ec22SJames Wright       PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, user, ceed_data, sgs_dd_setup_data));
6334c07ec22SJames Wright       break;
6344c07ec22SJames Wright   }
635ee1455b7SJames Wright 
6369c678832SJames Wright   // -- Create Operator to evalutate residual of SGS stress
637ad494f68SJames Wright   PetscCall(SgsSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data));
6389c678832SJames Wright 
639ad494f68SJames Wright   PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data));
640d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
64162b7942eSJames Wright }
642ee1455b7SJames Wright 
64342454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) {
644ee1455b7SJames Wright   PetscFunctionBeginUser;
645d949ddfcSJames Wright   if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS);
646b4c37c5cSJames Wright   Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed;
647ee1455b7SJames Wright 
648b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed));
649b87d60b3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed));
650ee1455b7SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx));
65141edf198SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx));
652b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx));
653b87d60b3SJames Wright   PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx));
654ee1455b7SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_sgs));
655b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs));
656b87d60b3SJames Wright   PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs));
657b87d60b3SJames 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));
658ee1455b7SJames Wright   PetscCall(PetscFree(sgs_dd_data));
659d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
660ee1455b7SJames Wright }
661