xref: /libCEED/examples/fluids/navierstokes.h (revision 07d5b98a8feba68a643190b8ea9bcdac5c3e6570)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 #pragma once
8 
9 #include <ceed-utils.h>
10 #include <ceed.h>
11 #include <mat-ceed.h>
12 #include <petscts.h>
13 #include <stdbool.h>
14 
15 #include "./include/petsc_ops.h"
16 #include "qfunctions/newtonian_types.h"
17 #include "qfunctions/stabilization_types.h"
18 
19 #if PETSC_VERSION_LT(3, 20, 0)
20 #error "PETSc v3.20 or later is required"
21 #endif
22 
23 #if PETSC_VERSION_LT(3, 21, 0)
24 #define DMSetCoordinateDisc(a, b, c) DMProjectCoordinates(a, b)
25 #define DMPlexFilter(a, b, c, d, e, f, g) DMPlexFilter(a, b, c, g)
26 #endif
27 
28 // -----------------------------------------------------------------------------
29 // Enums
30 // -----------------------------------------------------------------------------
31 
32 // Euler - test cases
33 typedef enum {
34   EULER_TEST_ISENTROPIC_VORTEX = 0,
35   EULER_TEST_1                 = 1,
36   EULER_TEST_2                 = 2,
37   EULER_TEST_3                 = 3,
38   EULER_TEST_4                 = 4,
39   EULER_TEST_5                 = 5,
40 } EulerTestType;
41 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1",      "test_2", "test_3", "test_4", "test_5",
42                                              "EulerTestType",     "EULER_TEST_", NULL};
43 
44 // Advection - Wind types
45 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL};
46 
47 // Advection - Initial Condition Types
48 static const char *const AdvectionICTypes[] = {"sphere", "cylinder", "cosine_hill", "skew", "AdvectionICType", "ADVECTIONIC_", NULL};
49 
50 // Advection - Bubble Continuity Types
51 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "cosine", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
52 
53 // Stabilization methods
54 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
55 
56 // Stabilization tau constants
57 static const char *const StabilizationTauTypes[] = {"Ctau", "AdvDiff_Shakib", "AdvDiff_Shakib_P", "StabilizationTauType", "STAB_TAU_", NULL};
58 
59 // Test mode type
60 typedef enum {
61   TESTTYPE_NONE           = 0,
62   TESTTYPE_SOLVER         = 1,
63   TESTTYPE_TURB_SPANSTATS = 2,
64   TESTTYPE_DIFF_FILTER    = 3,
65 } TestType;
66 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "diff_filter", "TestType", "TESTTYPE_", NULL};
67 
68 // Subgrid-Stress mode type
69 typedef enum {
70   SGS_MODEL_NONE        = 0,
71   SGS_MODEL_DATA_DRIVEN = 1,
72 } SGSModelType;
73 static const char *const SGSModelTypes[] = {"none", "data_driven", "SGSModelType", "SGS_MODEL_", NULL};
74 
75 // Mesh transformation type
76 typedef enum {
77   MESH_TRANSFORM_NONE      = 0,
78   MESH_TRANSFORM_PLATEMESH = 1,
79 } MeshTransformType;
80 static const char *const MeshTransformTypes[] = {"none", "platemesh", "MeshTransformType", "MESH_TRANSFORM_", NULL};
81 
82 static const char *const DifferentialFilterDampingFunctions[] = {
83     "none", "van_driest", "mms", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
84 
85 // -----------------------------------------------------------------------------
86 // Log Events
87 // -----------------------------------------------------------------------------
88 extern PetscLogEvent FLUIDS_CeedOperatorApply;
89 extern PetscLogEvent FLUIDS_CeedOperatorAssemble;
90 extern PetscLogEvent FLUIDS_CeedOperatorAssembleDiagonal;
91 extern PetscLogEvent FLUIDS_CeedOperatorAssemblePointBlockDiagonal;
92 extern PetscLogEvent FLUIDS_SmartRedis_Init;
93 extern PetscLogEvent FLUIDS_SmartRedis_Meta;
94 extern PetscLogEvent FLUIDS_SmartRedis_Train;
95 extern PetscLogEvent FLUIDS_TrainDataCompute;
96 extern PetscLogEvent FLUIDS_DifferentialFilter;
97 extern PetscLogEvent FLUIDS_VelocityGradientProjection;
98 PetscErrorCode       RegisterLogEvents();
99 
100 // -----------------------------------------------------------------------------
101 // Structs
102 // -----------------------------------------------------------------------------
103 // Structs declarations
104 typedef struct AppCtx_private   *AppCtx;
105 typedef struct CeedData_private *CeedData;
106 typedef struct User_private     *User;
107 typedef struct Units_private    *Units;
108 typedef struct SimpleBC_private *SimpleBC;
109 typedef struct Physics_private  *Physics;
110 
111 // Application context from user command line options
112 struct AppCtx_private {
113   // libCEED arguments
114   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
115   PetscInt degree;
116   PetscInt q_extra;
117   // Solver arguments
118   MatType amat_type;
119   // Post-processing arguments
120   PetscInt  checkpoint_interval;
121   PetscInt  viz_refine;
122   PetscInt  cont_steps;
123   PetscReal cont_time;
124   char      cont_file[PETSC_MAX_PATH_LEN];
125   char      cont_time_file[PETSC_MAX_PATH_LEN];
126   char      output_dir[PETSC_MAX_PATH_LEN];
127   PetscBool add_stepnum2bin;
128   PetscBool checkpoint_vtk;
129   // Problem type arguments
130   PetscFunctionList problems;
131   char              problem_name[PETSC_MAX_PATH_LEN];
132   // Test mode arguments
133   TestType    test_type;
134   PetscScalar test_tol;
135   char        test_file_path[PETSC_MAX_PATH_LEN];
136   // Turbulent spanwise statistics
137   PetscBool         turb_spanstats_enable;
138   PetscInt          turb_spanstats_collect_interval;
139   PetscInt          turb_spanstats_viewer_interval;
140   PetscViewer       turb_spanstats_viewer;
141   PetscViewerFormat turb_spanstats_viewer_format;
142   // Wall forces
143   struct {
144     PetscInt          num_wall;
145     PetscInt         *walls;
146     PetscViewer       viewer;
147     PetscViewerFormat viewer_format;
148     PetscBool         header_written;
149   } wall_forces;
150   // Subgrid Stress Model
151   SGSModelType sgs_model_type;
152   PetscBool    sgs_train_enable;
153   // Differential Filtering
154   PetscBool         diff_filter_monitor;
155   MeshTransformType mesh_transform_type;
156 };
157 
158 // libCEED data struct
159 struct CeedData_private {
160   CeedVector           x_coord, q_data;
161   CeedBasis            basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur;
162   CeedElemRestriction  elem_restr_x, elem_restr_q, elem_restr_qd_i;
163   CeedOperator         op_setup_vol;
164   OperatorApplyContext op_ics_ctx;
165   CeedQFunction        qf_setup_vol, qf_ics, qf_rhs_vol, qf_ifunction_vol, qf_setup_sur, qf_apply_inflow, qf_apply_inflow_jacobian, qf_apply_outflow,
166       qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian, qf_apply_slip, qf_apply_slip_jacobian;
167 };
168 
169 typedef struct {
170   DM                    dm;
171   PetscSF               sf;  // For communicating child data to parents
172   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
173   PetscInt              num_comp_stats;
174   Vec                   Child_Stats_loc, Parent_Stats_loc;
175   KSP                   ksp;         // For the L^2 projection solve
176   CeedScalar            span_width;  // spanwise width of the child domain
177   PetscBool             do_mms_test;
178   OperatorApplyContext  mms_error_ctx;
179   CeedContextFieldLabel solution_time_label, previous_time_label;
180 } SpanStatsData;
181 
182 typedef struct {
183   DM                   dm;
184   PetscInt             num_comp;
185   OperatorApplyContext l2_rhs_ctx;
186   KSP                  ksp;
187 } *NodalProjectionData;
188 
189 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
190 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
191 typedef struct {
192   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
193   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
194   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
195   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
196   SgsDDNodalStressEval      sgs_nodal_eval;
197   SgsDDNodalStressInference sgs_nodal_inference;
198   void                     *sgs_nodal_inference_ctx;
199   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
200 } *SgsDDData;
201 
202 typedef struct {
203   DM                   dm_dd_training;
204   PetscInt             num_comp_dd_inputs, write_data_interval;
205   OperatorApplyContext op_training_data_calc_ctx;
206   NodalProjectionData  filtered_grad_velo_proj;
207   size_t               training_data_array_dims[2];
208   PetscBool            overwrite_training_data;
209 } *SGS_DD_TrainingData;
210 
211 typedef struct {
212   DM                   dm_filter;
213   PetscInt             num_filtered_fields;
214   CeedInt             *num_field_components;
215   PetscInt             field_prim_state, field_velo_prod;
216   OperatorApplyContext op_rhs_ctx;
217   KSP                  ksp;
218   PetscBool            do_mms_test;
219 } *DiffFilterData;
220 
221 typedef struct {
222   void    *client;
223   char     rank_id_name[16];
224   PetscInt collocated_database_num_ranks;
225 } *SmartSimData;
226 
227 // PETSc user data
228 struct User_private {
229   MPI_Comm             comm;
230   DM                   dm;
231   DM                   dm_viz;
232   Mat                  interp_viz;
233   Ceed                 ceed;
234   Units                units;
235   Vec                  Q_loc, Q_dot_loc;
236   Physics              phys;
237   AppCtx               app_ctx;
238   CeedVector           q_ceed, q_dot_ceed, g_ceed, x_ceed;
239   CeedOperator         op_rhs_vol, op_ifunction_vol, op_ifunction;
240   Mat                  mat_ijacobian;
241   KSP                  mass_ksp;
242   OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx;
243   CeedScalar           time_bc_set;
244   SpanStatsData        spanstats;
245   NodalProjectionData  grad_velo_proj;
246   SgsDDData            sgs_dd_data;
247   DiffFilterData       diff_filter;
248   SmartSimData         smartsim;
249   SGS_DD_TrainingData  sgs_dd_train;
250 };
251 
252 // Units
253 struct Units_private {
254   // fundamental units
255   PetscScalar meter;
256   PetscScalar kilogram;
257   PetscScalar second;
258   PetscScalar Kelvin;
259   // derived units
260   PetscScalar Pascal;
261   PetscScalar J_per_kg_K;
262   PetscScalar m_per_squared_s;
263   PetscScalar W_per_m_K;
264   PetscScalar Joule;
265 };
266 
267 // Boundary conditions
268 struct SimpleBC_private {
269   PetscInt num_wall,  // Number of faces with wall BCs
270       wall_comps[5],  // An array of constrained component numbers
271       num_comps,
272       num_symmetry[3],  // Number of faces with symmetry BCs
273       num_inflow, num_outflow, num_freestream, num_slip;
274   PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16];
275 };
276 
277 // Struct that contains all enums and structs used for the physics of all problems
278 struct Physics_private {
279   PetscBool             implicit;
280   StateVariable         state_var;
281   CeedContextFieldLabel solution_time_label;
282   CeedContextFieldLabel stg_solution_time_label;
283   CeedContextFieldLabel timestep_size_label;
284   CeedContextFieldLabel ics_time_label;
285   CeedContextFieldLabel ijacobian_time_shift_label;
286 };
287 
288 typedef struct {
289   CeedQFunctionUser    qfunction;
290   const char          *qfunction_loc;
291   CeedQFunctionContext qfunction_context;
292 } ProblemQFunctionSpec;
293 
294 // Problem specific data
295 typedef struct ProblemData_private ProblemData;
296 struct ProblemData_private {
297   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
298   CeedScalar           dm_scale;
299   ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow,
300       apply_freestream, apply_slip, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
301   bool      non_zero_time;
302   PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian;
303   PetscErrorCode (*print_info)(User, ProblemData *, AppCtx);
304 };
305 
306 extern int FreeContextPetsc(void *);
307 
308 // -----------------------------------------------------------------------------
309 // Set up problems
310 // -----------------------------------------------------------------------------
311 // Set up function for each problem
312 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
313 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
314 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
315 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
316 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
317 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
318 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
319 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
320 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
321 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc);
322 
323 // Print function for each problem
324 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData *problem, AppCtx app_ctx);
325 
326 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData *problem, AppCtx app_ctx);
327 
328 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData *problem, AppCtx app_ctx);
329 
330 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData *problem, AppCtx app_ctx);
331 
332 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData *problem, AppCtx app_ctx);
333 
334 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData *problem, MPI_Comm comm);
335 
336 // -----------------------------------------------------------------------------
337 // libCEED functions
338 // -----------------------------------------------------------------------------
339 // Utility function to create local CEED restriction
340 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
341                                          CeedElemRestriction *elem_restr);
342 
343 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
344                                                CeedElemRestriction *restriction);
345 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
346                                                          CeedElemRestriction *restriction);
347 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
348                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
349 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
350                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
351 
352 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
353 
354 // Utility function to create CEED Composite Operator for the entire domain
355 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
356                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
357                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian);
358 
359 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc);
360 
361 // -----------------------------------------------------------------------------
362 // Time-stepping functions
363 // -----------------------------------------------------------------------------
364 // Create KSP to solve the inverse mass operator for explicit time stepping schemes
365 PetscErrorCode CreateKSPMassOperator(User user, CeedData ceed_data);
366 
367 // RHS (Explicit time-stepper) function setup
368 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
369 
370 // Implicit time-stepper function setup
371 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
372 
373 // User provided TS Monitor
374 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
375 
376 // TS: Create, setup, and solve
377 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts);
378 
379 // Update Boundary Values when time has changed
380 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
381 
382 // -----------------------------------------------------------------------------
383 // Setup DM
384 // -----------------------------------------------------------------------------
385 // Create mesh
386 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm);
387 
388 // Set up DM
389 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
390 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
391                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
392 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
393 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
394                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
395 
396 // Refine DM for high-order viz
397 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys);
398 
399 // -----------------------------------------------------------------------------
400 // Process command line options
401 // -----------------------------------------------------------------------------
402 // Register problems to be available on the command line
403 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
404 
405 // Process general command line options
406 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
407 
408 // -----------------------------------------------------------------------------
409 // Miscellaneous utility functions
410 // -----------------------------------------------------------------------------
411 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
412                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
413                                       CeedVector *inv_multiplicity);
414 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
415 
416 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
417                                                   Vec grad_FVM);
418 
419 // Compare reference solution values with current test run for CI
420 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
421 
422 // Get error for problems with exact solutions
423 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
424 
425 // Post-processing
426 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time);
427 
428 // -- Gather initial Q values in case of continuation of simulation
429 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
430 
431 // Record boundary values from initial condition
432 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
433 
434 // Versioning token for binary checkpoints
435 extern const PetscInt32 FLUIDS_FILE_TOKEN;  // for backwards compatibility
436 extern const PetscInt32 FLUIDS_FILE_TOKEN_32;
437 extern const PetscInt32 FLUIDS_FILE_TOKEN_64;
438 
439 // Create appropriate mass qfunction based on number of components N
440 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
441 
442 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
443 
444 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
445                                  FILE **fp);
446 
447 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
448 
449 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
450 
451 // -----------------------------------------------------------------------------
452 // Turbulence Statistics Collection Functions
453 // -----------------------------------------------------------------------------
454 
455 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
456 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
457 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data);
458 
459 // -----------------------------------------------------------------------------
460 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
461 // -----------------------------------------------------------------------------
462 
463 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
464 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
465 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc);
466 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, StateVariable state_var_input,
467                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
468 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
469 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
470                                                         CeedVector *grid_aniso_vector);
471 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
472                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
473 
474 // -----------------------------------------------------------------------------
475 // Boundary Condition Related Functions
476 // -----------------------------------------------------------------------------
477 
478 // Setup StrongBCs that use QFunctions
479 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc);
480 
481 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
482 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
483 PetscErrorCode SlipBCSetup(ProblemData *problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
484 
485 // -----------------------------------------------------------------------------
486 // Differential Filtering Functions
487 // -----------------------------------------------------------------------------
488 
489 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
490 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
491 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
492 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
493 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData *problem);
494 
495 // -----------------------------------------------------------------------------
496 // SGS Data-Driven Training via SmartSim
497 // -----------------------------------------------------------------------------
498 PetscErrorCode SmartSimSetup(User user);
499 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
500 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem);
501 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
502 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
503 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
504