xref: /honee/include/navierstokes.h (revision df29e1eeecb4505f1bf77a7dc8798babc49347ab)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 #pragma once
4 
5 #include <ceed.h>
6 #include <bc_definition.h>
7 #include <honee-file.h>
8 #include <log_events.h>
9 #include <mat-ceed.h>
10 #include <petsc-ceed-utils.h>
11 #include <petscts.h>
12 #include <stdbool.h>
13 
14 #include <petsc_ops.h>
15 #include "../qfunctions/newtonian_types.h"
16 
17 #if PETSC_VERSION_LT(3, 21, 0)
18 #error "PETSc v3.21 or later is required"
19 #endif
20 
21 // -----------------------------------------------------------------------------
22 // Enums
23 // -----------------------------------------------------------------------------
24 
25 // Euler - test cases
26 typedef enum {
27   EULER_TEST_ISENTROPIC_VORTEX = 0,
28   EULER_TEST_1                 = 1,
29   EULER_TEST_2                 = 2,
30   EULER_TEST_3                 = 3,
31   EULER_TEST_4                 = 4,
32   EULER_TEST_5                 = 5,
33 } EulerTestType;
34 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL};
35 
36 // Advection - Wind types
37 static const char *const AdvDifWindTypes[] = {"ROTATION", "TRANSLATION", "BOUNDARY_LAYER", "WindType", "ADVDIF_WIND_", NULL};
38 
39 // Advection - Initial Condition Types
40 static const char *const AdvDifICTypes[] = {"SPHERE",         "CYLINDER",        "COSINE_HILL", "SKEW", "WAVE",
41                                             "BOUNDARY_LAYER", "AdvectionICType", "ADVDIF_IC_",  NULL};
42 
43 // Advection - Wave Types
44 static const char *const AdvDifWaveTypes[] = {"SINE", "SQUARE", "AdvDifWaveType", "ADVDIF_WAVE_", NULL};
45 
46 // Advection - Bubble Continuity Types
47 static const char *const AdvDifBubbleContinuityTypes[] = {
48     "SMOOTH", "BACK_SHARP", "THICK", "COSINE", "BubbleContinuityType", "ADVDIF_BUBBLE_CONTINUITY_", NULL};
49 
50 // Stabilization methods
51 static const char *const StabilizationTypes[] = {"NONE", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
52 
53 // Stabilization tau constants
54 static const char *const StabilizationTauTypes[] = {"CTAU", "ADVDIFF_SHAKIB", "StabilizationTauType", "STAB_TAU_", NULL};
55 
56 // Test mode type
57 typedef enum {
58   TESTTYPE_NONE           = 0,
59   TESTTYPE_SOLVER         = 1,
60   TESTTYPE_TURB_SPANSTATS = 2,
61   TESTTYPE_DIFF_FILTER    = 3,
62 } TestType;
63 static const char *const TestTypes[] = {"NONE", "SOLVER", "TURB_SPANSTATS", "DIFF_FILTER", "TestType", "TESTTYPE_", NULL};
64 
65 // Subgrid-Stress mode type
66 typedef enum {
67   SGS_MODEL_NONE        = 0,
68   SGS_MODEL_DATA_DRIVEN = 1,
69 } SGSModelType;
70 static const char *const SGSModelTypes[] = {"NONE", "DATA_DRIVEN", "SGSModelType", "SGS_MODEL_", NULL};
71 
72 // Subgrid-Stress mode type
73 typedef enum {
74   SGS_MODEL_DD_FUSED           = 0,
75   SGS_MODEL_DD_SEQENTIAL_CEED  = 1,
76   SGS_MODEL_DD_SEQENTIAL_TORCH = 2,
77 } SGSModelDDImplementation;
78 static const char *const SGSModelDDImplementations[] = {"FUSED", "SEQUENTIAL_CEED", "SEQUENTIAL_TORCH", "SGSModelDDImplementation", "SGS_MODEL_DD_",
79                                                         NULL};
80 
81 // Mesh transformation type
82 typedef enum {
83   MESH_TRANSFORM_NONE      = 0,
84   MESH_TRANSFORM_PLATEMESH = 1,
85 } MeshTransformType;
86 static const char *const MeshTransformTypes[] = {"NONE", "PLATEMESH", "MeshTransformType", "MESH_TRANSFORM_", NULL};
87 
88 static const char *const DifferentialFilterDampingFunctions[] = {
89     "NONE", "VAN_DRIEST", "MMS", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
90 
91 static const char *const DivDiffFluxProjectionMethods[] = {"NONE", "DIRECT", "INDIRECT", "DivDiffFluxProjectionMethod", "DIV_DIFF_FLUX_PROJ_", NULL};
92 
93 // -----------------------------------------------------------------------------
94 // Structs
95 // -----------------------------------------------------------------------------
96 // Structs declarations
97 typedef struct AppCtx_private      *AppCtx;
98 typedef struct CeedData_private    *CeedData;
99 typedef struct Honee_private       *Honee;
100 typedef struct Units_private       *Units;
101 typedef struct SimpleBC_private    *SimpleBC;
102 typedef struct Physics_private     *Physics;
103 typedef struct ProblemData_private *ProblemData;
104 
105 // Application context from user command line options
106 struct AppCtx_private {
107   // libCEED arguments
108   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
109   PetscInt degree;
110   PetscInt q_extra;
111   // Solver arguments
112   MatType amat_type;
113   // Post-processing arguments
114   PetscInt  checkpoint_interval;
115   PetscInt  viz_refine;
116   PetscInt  cont_steps;
117   PetscReal cont_time;
118   char      cont_file[PETSC_MAX_PATH_LEN];
119   char      cont_time_file[PETSC_MAX_PATH_LEN];
120   char      output_dir[PETSC_MAX_PATH_LEN];
121   PetscBool add_stepnum2bin;
122   PetscBool checkpoint_vtk;
123   // Problem type arguments
124   PetscFunctionList problems;
125   char              problem_name[PETSC_MAX_PATH_LEN];
126   // Test mode arguments
127   TestType    test_type;
128   PetscScalar test_tol;
129   char        test_file_path[PETSC_MAX_PATH_LEN];
130   // Turbulent spanwise statistics
131   PetscBool         turb_spanstats_enable;
132   PetscInt          turb_spanstats_collect_interval;
133   PetscInt          turb_spanstats_viewer_interval;
134   PetscViewer       turb_spanstats_viewer;
135   PetscViewerFormat turb_spanstats_viewer_format;
136   // Wall forces
137   struct {
138     PetscInt          num_wall;
139     PetscInt         *walls;
140     PetscViewer       viewer;
141     PetscViewerFormat viewer_format;
142     PetscBool         header_written;
143   } wall_forces;
144   // Subgrid Stress Model
145   SGSModelType sgs_model_type;
146   PetscBool    sgs_train_enable;
147   // Differential Filtering
148   PetscBool         diff_filter_monitor;
149   MeshTransformType mesh_transform_type;
150   // Divergence of Diffusive Flux Projection
151   DivDiffFluxProjectionMethod divFdiffproj_method;
152 };
153 
154 // libCEED data struct
155 struct CeedData_private {
156   CeedVector           x_coord;
157   CeedBasis            basis_x, basis_q;
158   CeedElemRestriction  elem_restr_x, elem_restr_q;
159   OperatorApplyContext op_ics_ctx;
160 };
161 
162 typedef struct {
163   DM                    dm;
164   PetscSF               sf;  // For communicating child data to parents
165   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
166   PetscInt              num_comp_stats;
167   Vec                   Child_Stats_loc, Parent_Stats_loc;
168   KSP                   ksp;         // For the L^2 projection solve
169   CeedScalar            span_width;  // spanwise width of the child domain
170   PetscBool             do_mms_test;
171   OperatorApplyContext  mms_error_ctx;
172   CeedContextFieldLabel solution_time_label, previous_time_label;
173 } SpanStatsData;
174 
175 typedef struct {
176   DM                   dm;
177   PetscInt             num_comp;
178   OperatorApplyContext l2_rhs_ctx;
179   KSP                  ksp;
180 } *NodalProjectionData;
181 
182 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData;
183 
184 struct DivDiffFluxProjectionData_private {
185   PetscInt                    num_diff_flux_comps;
186   DivDiffFluxProjectionMethod method;
187   NodalProjectionData         projection;
188 
189   // CeedOperator Objects
190   CeedElemRestriction elem_restr_div_diff_flux;
191   CeedBasis           basis_div_diff_flux;
192   CeedEvalMode        eval_mode_div_diff_flux;
193   CeedVector          div_diff_flux_ceed;
194 
195   // Problem specific setup functions
196   PetscErrorCode (*CreateRHSOperator_Direct)(Honee, CeedData, DivDiffFluxProjectionData, CeedOperator *);
197   PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, CeedData, DivDiffFluxProjectionData, CeedOperator *);
198 
199   // Only used for direct method:
200   Vec          DivDiffFlux_loc;
201   PetscMemType DivDiffFlux_memtype;
202   PetscBool    ceed_vec_has_array;
203 
204   // Only used for indirect method:
205   OperatorApplyContext calc_div_diff_flux;
206 };
207 
208 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
209 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
210 typedef struct {
211   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
212   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
213   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
214   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
215   SgsDDNodalStressEval      sgs_nodal_eval;
216   SgsDDNodalStressInference sgs_nodal_inference;
217   void                     *sgs_nodal_inference_ctx;
218   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
219 } *SgsDDData;
220 
221 typedef struct {
222   DM                   dm_dd_training;
223   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
224   PetscScalar          filter_widths[16];
225   OperatorApplyContext op_training_data_calc_ctx;
226   NodalProjectionData  filtered_grad_velo_proj;
227   size_t               training_data_array_dims[2];
228   PetscBool            overwrite_training_data;
229 } *SGS_DD_TrainingData;
230 
231 typedef struct {
232   DM                    dm_filter;
233   PetscInt              num_filtered_fields;
234   CeedInt              *num_field_components;
235   PetscInt              field_prim_state, field_velo_prod;
236   OperatorApplyContext  op_rhs_ctx;
237   KSP                   ksp;
238   PetscObjectState      X_loc_state;
239   PetscBool             do_mms_test;
240   CeedContextFieldLabel filter_width_scaling_label;
241 } *DiffFilterData;
242 
243 typedef struct {
244   void    *client;
245   char     rank_id_name[16];
246   PetscInt collocated_database_num_ranks;
247 } *SmartSimData;
248 
249 // PETSc user data
250 struct Honee_private {
251   MPI_Comm                  comm;
252   DM                        dm;
253   DM                        dm_viz;
254   Mat                       interp_viz;
255   Ceed                      ceed;
256   Units                     units;
257   Vec                       Q_loc, Q_dot_loc;
258   Physics                   phys;
259   AppCtx                    app_ctx;
260   CeedVector                q_ceed, q_dot_ceed, g_ceed, x_ceed;
261   CeedOperator              op_ifunction;
262   Mat                       mat_ijacobian;
263   KSP                       mass_ksp;
264   OperatorApplyContext      op_rhs_ctx, op_strong_bc_ctx;
265   CeedScalar                time_bc_set;
266   SpanStatsData             spanstats;
267   NodalProjectionData       grad_velo_proj;
268   SgsDDData                 sgs_dd_data;
269   DiffFilterData            diff_filter;
270   SmartSimData              smartsim;
271   SGS_DD_TrainingData       sgs_dd_train;
272   DivDiffFluxProjectionData diff_flux_proj;
273 };
274 
275 // Units
276 struct Units_private {
277   // fundamental units
278   PetscScalar meter;
279   PetscScalar kilogram;
280   PetscScalar second;
281   PetscScalar Kelvin;
282   // derived units
283   PetscScalar Pascal;
284   PetscScalar J_per_kg_K;
285   PetscScalar m_per_squared_s;
286   PetscScalar W_per_m_K;
287   PetscScalar Joule;
288 };
289 
290 // Boundary conditions
291 struct SimpleBC_private {
292   PetscInt num_inflow, num_outflow, num_freestream, num_slip;
293   PetscInt inflows[16], outflows[16], freestreams[16], slips[16];
294 };
295 
296 // Struct that contains all enums and structs used for the physics of all problems
297 struct Physics_private {
298   PetscBool             implicit;
299   StateVariable         state_var;
300   CeedContextFieldLabel solution_time_label;
301   CeedContextFieldLabel stg_solution_time_label;
302   CeedContextFieldLabel timestep_size_label;
303   CeedContextFieldLabel ics_time_label;
304 };
305 
306 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
307 
308 typedef struct {
309   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
310   const char          *qf_loc;       // !< Absolute path to QFunction source file
311   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
312 } ProblemQFunctionSpec;
313 
314 // Problem specific data
315 struct ProblemData_private {
316   CeedInt              jac_data_size_vol, jac_data_size_sur;
317   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
318       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
319   bool          compute_exact_solution_error;
320   PetscBool     set_bc_from_ics, use_strong_bc_ceed;
321   PetscCount    num_bc_defs;
322   BCDefinition *bc_defs;
323   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
324   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
325 };
326 
327 extern int FreeContextPetsc(void *);
328 
329 // -----------------------------------------------------------------------------
330 // Set up problems
331 // -----------------------------------------------------------------------------
332 // Set up function for each problem
333 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
334 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
335 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
336 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
337 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
338 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
339 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
340 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
341 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
342 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
343 
344 // Print function for each problem
345 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx);
346 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx);
347 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx);
348 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx);
349 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx);
350 
351 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
352 
353 // -----------------------------------------------------------------------------
354 // libCEED functions
355 // -----------------------------------------------------------------------------
356 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
357                                                CeedElemRestriction *restriction);
358 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
359                                                          CeedElemRestriction *restriction);
360 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
361                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
362 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
363                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
364 
365 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
366 PetscErrorCode CreateCoordinateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis);
367 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face,
368                                                          CeedBasis *basis);
369 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field,
370                                                CeedBasis *basis);
371 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name);
372 PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array);
373 
374 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
375 
376 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
377                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
378 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
379 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
380                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
381 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
382 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
383 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
384                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
385 PetscErrorCode QDataClearStoredData();
386 // -----------------------------------------------------------------------------
387 // Time-stepping functions
388 // -----------------------------------------------------------------------------
389 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
390 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
391 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
392 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts);
393 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
394 
395 // -----------------------------------------------------------------------------
396 // Setup DM
397 // -----------------------------------------------------------------------------
398 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
399 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
400 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys);
401 
402 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
403                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
404 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
405 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
406                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
407 
408 // -----------------------------------------------------------------------------
409 // Process command line options
410 // -----------------------------------------------------------------------------
411 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
412 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
413 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
414 
415 // -----------------------------------------------------------------------------
416 // Miscellaneous utility functions
417 // -----------------------------------------------------------------------------
418 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
419                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
420                                       CeedVector *inv_multiplicity);
421 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
422 
423 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
424                                                   Vec grad_FVM);
425 
426 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
427 PetscErrorCode PrintError(CeedData ceed_data, DM dm, Honee honee, Vec Q, PetscScalar final_time);
428 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
429 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
430 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
431 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
432 
433 // -----------------------------------------------------------------------------
434 // Turbulence Statistics Collection Functions
435 // -----------------------------------------------------------------------------
436 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem);
437 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
438 PetscErrorCode TurbulenceStatisticsDestroy(Honee honee, CeedData ceed_data);
439 
440 // -----------------------------------------------------------------------------
441 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
442 // -----------------------------------------------------------------------------
443 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem);
444 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
445 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
446 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem, StateVariable state_var_input,
447                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
448 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
449 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
450                                                         CeedVector *grid_aniso_vector);
451 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
452                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
453 
454 // -----------------------------------------------------------------------------
455 // Boundary Condition Related Functions
456 // -----------------------------------------------------------------------------
457 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, Honee honee, ProblemData problem, SimpleBC bc);
458 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
459 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
460 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
461 
462 // -----------------------------------------------------------------------------
463 // Differential Filtering Functions
464 // -----------------------------------------------------------------------------
465 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem);
466 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
467 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
468 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
469 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
470 
471 // -----------------------------------------------------------------------------
472 // SGS Data-Driven Training via SmartSim
473 // -----------------------------------------------------------------------------
474 PetscErrorCode SmartSimSetup(Honee honee);
475 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
476 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem);
477 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
478 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
479 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
480 
481 // -----------------------------------------------------------------------------
482 // Divergence of Diffusive Flux Projection
483 // -----------------------------------------------------------------------------
484 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj);
485 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
486                                                          CeedVector *vector, CeedEvalMode *eval_mode);
487 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, CeedData ceed_data, DivDiffFluxProjectionData diff_flux_proj);
488 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
489 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
490