xref: /honee/include/navierstokes.h (revision 481d14cbe61cedfaf9e39a526b7a6439011a8661)
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, 22, 0)
18 #error "PETSc v3.22 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 Honee_private       *Honee;
99 typedef struct Units_private       *Units;
100 typedef struct SimpleBC_private    *SimpleBC;
101 typedef struct Physics_private     *Physics;
102 typedef struct ProblemData_private *ProblemData;
103 
104 // Application context from user command line options
105 struct AppCtx_private {
106   // libCEED arguments
107   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
108   PetscInt degree;
109   PetscInt q_extra;
110   // Solver arguments
111   MatType amat_type;
112   // Post-processing arguments
113   PetscInt  checkpoint_interval;
114   PetscInt  viz_refine;
115   PetscBool use_continue_file;
116   PetscInt  cont_steps;
117   PetscReal cont_time;
118   char      cont_file[PETSC_MAX_PATH_LEN];
119   char      output_dir[PETSC_MAX_PATH_LEN];
120   PetscBool add_stepnum2bin;
121   PetscBool checkpoint_vtk;
122   // Problem type arguments
123   PetscFunctionList problems;
124   char              problem_name[PETSC_MAX_PATH_LEN];
125   // Test mode arguments
126   TestType    test_type;
127   PetscScalar test_tol;
128   char        test_file_path[PETSC_MAX_PATH_LEN];
129   // Turbulent spanwise statistics
130   PetscBool         turb_spanstats_enable;
131   PetscInt          turb_spanstats_collect_interval;
132   PetscInt          turb_spanstats_viewer_interval;
133   PetscViewer       turb_spanstats_viewer;
134   PetscViewerFormat turb_spanstats_viewer_format;
135   // Wall forces
136   struct {
137     PetscInt          num_wall;
138     PetscInt         *walls;
139     PetscViewer       viewer;
140     PetscViewerFormat viewer_format;
141     PetscBool         header_written;
142   } wall_forces;
143   // Subgrid Stress Model
144   SGSModelType sgs_model_type;
145   PetscBool    sgs_train_enable;
146   // Differential Filtering
147   PetscBool         diff_filter_monitor;
148   MeshTransformType mesh_transform_type;
149   // Divergence of Diffusive Flux Projection
150   DivDiffFluxProjectionMethod divFdiffproj_method;
151 };
152 
153 typedef struct {
154   DM                    dm;
155   PetscSF               sf;  // For communicating child data to parents
156   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
157   PetscInt              num_comp_stats;
158   Vec                   Child_Stats_loc, Parent_Stats_loc;
159   KSP                   ksp;         // For the L^2 projection solve
160   CeedScalar            span_width;  // spanwise width of the child domain
161   PetscBool             do_mms_test;
162   OperatorApplyContext  mms_error_ctx;
163   CeedContextFieldLabel solution_time_label, previous_time_label;
164 } SpanStatsData;
165 
166 typedef struct {
167   DM                   dm;
168   PetscInt             num_comp;
169   OperatorApplyContext l2_rhs_ctx;
170   KSP                  ksp;
171 } *NodalProjectionData;
172 
173 typedef struct DivDiffFluxProjectionData_private *DivDiffFluxProjectionData;
174 
175 struct DivDiffFluxProjectionData_private {
176   PetscInt                    num_diff_flux_comps;
177   DivDiffFluxProjectionMethod method;
178   NodalProjectionData         projection;
179 
180   // CeedOperator Objects
181   CeedElemRestriction elem_restr_div_diff_flux;
182   CeedBasis           basis_div_diff_flux;
183   CeedEvalMode        eval_mode_div_diff_flux;
184   CeedVector          div_diff_flux_ceed;
185 
186   // Problem specific setup functions
187   PetscErrorCode (*CreateRHSOperator_Direct)(Honee, DivDiffFluxProjectionData, CeedOperator *);
188   PetscErrorCode (*CreateRHSOperator_Indirect)(Honee, DivDiffFluxProjectionData, CeedOperator *);
189 
190   // Only used for direct method:
191   Vec          DivDiffFlux_loc;
192   PetscMemType DivDiffFlux_memtype;
193   PetscBool    ceed_vec_has_array;
194 
195   // Only used for indirect method:
196   OperatorApplyContext calc_div_diff_flux;
197 };
198 
199 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
200 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
201 typedef struct {
202   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
203   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
204   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
205   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
206   SgsDDNodalStressEval      sgs_nodal_eval;
207   SgsDDNodalStressInference sgs_nodal_inference;
208   void                     *sgs_nodal_inference_ctx;
209   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
210 } *SgsDDData;
211 
212 typedef struct {
213   DM                   dm_dd_training;
214   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
215   PetscScalar          filter_widths[16];
216   OperatorApplyContext op_training_data_calc_ctx;
217   NodalProjectionData  filtered_grad_velo_proj;
218   size_t               training_data_array_dims[2];
219   PetscBool            overwrite_training_data;
220 } *SGS_DD_TrainingData;
221 
222 typedef struct {
223   DM                    dm_filter;
224   PetscInt              num_filtered_fields;
225   CeedInt              *num_field_components;
226   PetscInt              field_prim_state, field_velo_prod;
227   OperatorApplyContext  op_rhs_ctx;
228   KSP                   ksp;
229   PetscObjectState      X_loc_state;
230   PetscBool             do_mms_test;
231   CeedContextFieldLabel filter_width_scaling_label;
232 } *DiffFilterData;
233 
234 typedef struct {
235   void    *client;
236   char     rank_id_name[16];
237   PetscInt collocated_database_num_ranks;
238 } *SmartSimData;
239 
240 // PETSc user data
241 struct Honee_private {
242   MPI_Comm                  comm;
243   DM                        dm;
244   DM                        dm_viz;
245   Mat                       interp_viz;
246   Ceed                      ceed;
247   Units                     units;
248   Vec                       Q_loc, Q_dot_loc;
249   Physics                   phys;
250   AppCtx                    app_ctx;
251   CeedVector                q_ceed, q_dot_ceed, g_ceed, x_ceed;
252   CeedOperator              op_ifunction;
253   Mat                       mat_ijacobian;
254   KSP                       mass_ksp;
255   OperatorApplyContext      op_rhs_ctx, op_strong_bc_ctx;
256   CeedScalar                time_bc_set;
257   SpanStatsData             spanstats;
258   NodalProjectionData       grad_velo_proj;
259   SgsDDData                 sgs_dd_data;
260   DiffFilterData            diff_filter;
261   SmartSimData              smartsim;
262   SGS_DD_TrainingData       sgs_dd_train;
263   DivDiffFluxProjectionData diff_flux_proj;
264 
265   // old CeedData
266   CeedVector           x_coord;
267   CeedBasis            basis_x, basis_q;
268   CeedElemRestriction  elem_restr_x, elem_restr_q;
269   OperatorApplyContext op_ics_ctx;
270 };
271 
272 // Units
273 struct Units_private {
274   // fundamental units
275   PetscScalar meter;
276   PetscScalar kilogram;
277   PetscScalar second;
278   PetscScalar Kelvin;
279   // derived units
280   PetscScalar Pascal;
281   PetscScalar J_per_kg_K;
282   PetscScalar m_per_squared_s;
283   PetscScalar W_per_m_K;
284   PetscScalar Joule;
285 };
286 
287 // Boundary conditions
288 struct SimpleBC_private {
289   PetscInt num_inflow, num_outflow, num_freestream, num_slip;
290   PetscInt inflows[16], outflows[16], freestreams[16], slips[16];
291 };
292 
293 // Struct that contains all enums and structs used for the physics of all problems
294 struct Physics_private {
295   PetscBool             implicit;
296   StateVariable         state_var;
297   CeedContextFieldLabel solution_time_label;
298   CeedContextFieldLabel stg_solution_time_label;
299   CeedContextFieldLabel timestep_size_label;
300   CeedContextFieldLabel ics_time_label;
301 };
302 
303 PetscErrorCode BoundaryConditionSetUp(Honee honee, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
304 
305 typedef struct {
306   CeedQFunctionUser    qf_func_ptr;  // !< QFunction function pointer
307   const char          *qf_loc;       // !< Absolute path to QFunction source file
308   CeedQFunctionContext qfctx;        // !< QFunctionContext to attach to QFunction
309 } ProblemQFunctionSpec;
310 
311 // Problem specific data
312 struct ProblemData_private {
313   CeedInt              jac_data_size_vol, jac_data_size_sur;
314   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
315       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
316   bool          compute_exact_solution_error;
317   PetscBool     set_bc_from_ics, use_strong_bc_ceed;
318   PetscCount    num_bc_defs;
319   BCDefinition *bc_defs;
320   PetscErrorCode (*print_info)(Honee, ProblemData, AppCtx);
321   PetscErrorCode (*create_mass_operator)(Honee, CeedOperator *);
322 };
323 
324 extern int FreeContextPetsc(void *);
325 
326 // -----------------------------------------------------------------------------
327 // Set up problems
328 // -----------------------------------------------------------------------------
329 // Set up function for each problem
330 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
331 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
332 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
333 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
334 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
335 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
336 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
337 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
338 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
339 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
340 
341 // Print function for each problem
342 extern PetscErrorCode PRINT_NEWTONIAN(Honee honee, ProblemData problem, AppCtx app_ctx);
343 extern PetscErrorCode PRINT_EULER_VORTEX(Honee honee, ProblemData problem, AppCtx app_ctx);
344 extern PetscErrorCode PRINT_SHOCKTUBE(Honee honee, ProblemData problem, AppCtx app_ctx);
345 extern PetscErrorCode PRINT_ADVECTION(Honee honee, ProblemData problem, AppCtx app_ctx);
346 extern PetscErrorCode PRINT_ADVECTION2D(Honee honee, ProblemData problem, AppCtx app_ctx);
347 
348 PetscErrorCode PrintRunInfo(Honee honee, Physics phys_ctx, ProblemData problem, TS ts);
349 
350 // -----------------------------------------------------------------------------
351 // libCEED functions
352 // -----------------------------------------------------------------------------
353 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
354                                                CeedElemRestriction *restriction);
355 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
356                                                          CeedElemRestriction *restriction);
357 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
358                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
359 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
360                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
361 
362 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
363 PetscErrorCode CreateCoordinateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, CeedBasis *basis);
364 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face,
365                                                          CeedBasis *basis);
366 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field,
367                                                CeedBasis *basis);
368 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, char **face_label_name);
369 PetscErrorCode DMLabelCreateGlobalValueArray(DM dm, DMLabel label, PetscInt *num_values, PetscInt **value_array);
370 
371 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
372 
373 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
374                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
375 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
376 PetscErrorCode QDataBoundaryGet(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 QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
379 PetscErrorCode QDataBoundaryGradientGetNumComponents(DM dm, CeedInt *q_data_size);
380 PetscErrorCode QDataBoundaryGradientGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedVector x_coord,
381                                         CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
382 PetscErrorCode QDataClearStoredData();
383 // -----------------------------------------------------------------------------
384 // Time-stepping functions
385 // -----------------------------------------------------------------------------
386 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
387 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
388 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
389 PetscErrorCode TSSolve_NS(DM dm, Honee honee, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts);
390 PetscErrorCode UpdateBoundaryValues(Honee honee, Vec Q_loc, PetscReal t);
391 
392 // -----------------------------------------------------------------------------
393 // Setup DM
394 // -----------------------------------------------------------------------------
395 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
396 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
397 PetscErrorCode VizRefineDM(DM dm, Honee honee, ProblemData problem, SimpleBC bc, Physics phys);
398 
399 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
400                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
401 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
402 PetscErrorCode DMSetupByOrder_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 
405 // -----------------------------------------------------------------------------
406 // Process command line options
407 // -----------------------------------------------------------------------------
408 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
409 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
410 PetscErrorCode HoneeOptionsSetValueDefault(PetscOptions options, const char name[], const char value[]);
411 
412 // -----------------------------------------------------------------------------
413 // Miscellaneous utility functions
414 // -----------------------------------------------------------------------------
415 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
416                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
417                                       CeedVector *inv_multiplicity);
418 PetscErrorCode ICs_FixMultiplicity(DM dm, Honee honee, Vec Q_loc, Vec Q, CeedScalar time);
419 
420 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
421                                                   Vec grad_FVM);
422 
423 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
424 PetscErrorCode PrintError(DM dm, Honee honee, Vec Q, PetscScalar final_time);
425 PetscErrorCode PostProcess(TS ts, DM dm, ProblemData problem, Honee honee, Vec Q, PetscScalar final_time);
426 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
427 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
428 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
429 
430 // -----------------------------------------------------------------------------
431 // Turbulence Statistics Collection Functions
432 // -----------------------------------------------------------------------------
433 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, Honee honee, ProblemData problem);
434 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
435 PetscErrorCode TurbulenceStatisticsDestroy(Honee honee);
436 
437 // -----------------------------------------------------------------------------
438 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
439 // -----------------------------------------------------------------------------
440 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem);
441 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
442 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc);
443 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, Honee honee, ProblemData problem, StateVariable state_var_input,
444                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
445 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
446 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
447                                                         CeedVector *grid_aniso_vector);
448 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, Honee honee, CeedElemRestriction *elem_restr_grid_aniso,
449                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
450 
451 // -----------------------------------------------------------------------------
452 // Boundary Condition Related Functions
453 // -----------------------------------------------------------------------------
454 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem, SimpleBC bc);
455 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
456 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
457 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
458 
459 // -----------------------------------------------------------------------------
460 // Differential Filtering Functions
461 // -----------------------------------------------------------------------------
462 PetscErrorCode DifferentialFilterSetup(Ceed ceed, Honee honee, ProblemData problem);
463 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
464 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
465 PetscErrorCode DifferentialFilterApply(Honee honee, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
466 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
467 
468 // -----------------------------------------------------------------------------
469 // SGS Data-Driven Training via SmartSim
470 // -----------------------------------------------------------------------------
471 PetscErrorCode SmartSimSetup(Honee honee);
472 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
473 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem);
474 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
475 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
476 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
477 
478 // -----------------------------------------------------------------------------
479 // Divergence of Diffusive Flux Projection
480 // -----------------------------------------------------------------------------
481 PetscErrorCode DivDiffFluxProjectionCreate(Honee honee, PetscInt num_diff_flux_comps, DivDiffFluxProjectionData *pdiff_flux_proj);
482 PetscErrorCode DivDiffFluxProjectionGetOperatorFieldData(DivDiffFluxProjectionData diff_flux_proj, CeedElemRestriction *elem_restr, CeedBasis *basis,
483                                                          CeedVector *vector, CeedEvalMode *eval_mode);
484 PetscErrorCode DivDiffFluxProjectionSetup(Honee honee, DivDiffFluxProjectionData diff_flux_proj);
485 PetscErrorCode DivDiffFluxProjectionApply(DivDiffFluxProjectionData diff_flux_proj, Vec Q_loc);
486 PetscErrorCode DivDiffFluxProjectionDataDestroy(DivDiffFluxProjectionData diff_flux_proj);
487