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