xref: /honee/include/navierstokes.h (revision 21ba7ba458ca7fbd5416be1f0b0586f2885778ca)
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.h>
10 #include <bc_definition.h>
11 #include <log_events.h>
12 #include <mat-ceed.h>
13 #include <petsc-ceed-utils.h>
14 #include <petscts.h>
15 #include <stdbool.h>
16 
17 #include <petsc_ops.h>
18 #include "../qfunctions/newtonian_types.h"
19 
20 #if PETSC_VERSION_LT(3, 21, 0)
21 #error "PETSc v3.21 or later is required"
22 #endif
23 
24 // -----------------------------------------------------------------------------
25 // Enums
26 // -----------------------------------------------------------------------------
27 
28 // Euler - test cases
29 typedef enum {
30   EULER_TEST_ISENTROPIC_VORTEX = 0,
31   EULER_TEST_1                 = 1,
32   EULER_TEST_2                 = 2,
33   EULER_TEST_3                 = 3,
34   EULER_TEST_4                 = 4,
35   EULER_TEST_5                 = 5,
36 } EulerTestType;
37 static const char *const EulerTestTypes[] = {"ISENTROPIC_VORTEX", "1", "2", "3", "4", "5", "EulerTestType", "EULER_TEST_", NULL};
38 
39 // Advection - Wind types
40 static const char *const WindTypes[] = {"ROTATION", "TRANSLATION", "WindType", "WIND_", NULL};
41 
42 // Advection - Initial Condition Types
43 static const char *const AdvectionICTypes[] = {"SPHERE", "CYLINDER", "COSINE_HILL", "SKEW", "AdvectionICType", "ADVECTIONIC_", 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", "ADVDIFF_SHAKIB_P", "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 // -----------------------------------------------------------------------------
90 // Structs
91 // -----------------------------------------------------------------------------
92 // Structs declarations
93 typedef struct AppCtx_private      *AppCtx;
94 typedef struct CeedData_private    *CeedData;
95 typedef struct User_private        *User;
96 typedef struct Units_private       *Units;
97 typedef struct SimpleBC_private    *SimpleBC;
98 typedef struct Physics_private     *Physics;
99 typedef struct ProblemData_private *ProblemData;
100 
101 // Application context from user command line options
102 struct AppCtx_private {
103   // libCEED arguments
104   char     ceed_resource[PETSC_MAX_PATH_LEN];  // libCEED backend
105   PetscInt degree;
106   PetscInt q_extra;
107   // Solver arguments
108   MatType amat_type;
109   // Post-processing arguments
110   PetscInt  checkpoint_interval;
111   PetscInt  viz_refine;
112   PetscInt  cont_steps;
113   PetscReal cont_time;
114   char      cont_file[PETSC_MAX_PATH_LEN];
115   char      cont_time_file[PETSC_MAX_PATH_LEN];
116   char      output_dir[PETSC_MAX_PATH_LEN];
117   PetscBool add_stepnum2bin;
118   PetscBool checkpoint_vtk;
119   // Problem type arguments
120   PetscFunctionList problems;
121   char              problem_name[PETSC_MAX_PATH_LEN];
122   // Test mode arguments
123   TestType    test_type;
124   PetscScalar test_tol;
125   char        test_file_path[PETSC_MAX_PATH_LEN];
126   // Turbulent spanwise statistics
127   PetscBool         turb_spanstats_enable;
128   PetscInt          turb_spanstats_collect_interval;
129   PetscInt          turb_spanstats_viewer_interval;
130   PetscViewer       turb_spanstats_viewer;
131   PetscViewerFormat turb_spanstats_viewer_format;
132   // Wall forces
133   struct {
134     PetscInt          num_wall;
135     PetscInt         *walls;
136     PetscViewer       viewer;
137     PetscViewerFormat viewer_format;
138     PetscBool         header_written;
139   } wall_forces;
140   // Subgrid Stress Model
141   SGSModelType sgs_model_type;
142   PetscBool    sgs_train_enable;
143   // Differential Filtering
144   PetscBool         diff_filter_monitor;
145   MeshTransformType mesh_transform_type;
146 };
147 
148 // libCEED data struct
149 struct CeedData_private {
150   CeedVector           x_coord, q_data;
151   CeedBasis            basis_x, basis_q;
152   CeedElemRestriction  elem_restr_x, elem_restr_q, elem_restr_qd_i;
153   OperatorApplyContext op_ics_ctx;
154 };
155 
156 typedef struct {
157   DM                    dm;
158   PetscSF               sf;  // For communicating child data to parents
159   OperatorApplyContext  op_stats_collect_ctx, op_proj_rhs_ctx;
160   PetscInt              num_comp_stats;
161   Vec                   Child_Stats_loc, Parent_Stats_loc;
162   KSP                   ksp;         // For the L^2 projection solve
163   CeedScalar            span_width;  // spanwise width of the child domain
164   PetscBool             do_mms_test;
165   OperatorApplyContext  mms_error_ctx;
166   CeedContextFieldLabel solution_time_label, previous_time_label;
167 } SpanStatsData;
168 
169 typedef struct {
170   DM                   dm;
171   PetscInt             num_comp;
172   OperatorApplyContext l2_rhs_ctx;
173   KSP                  ksp;
174 } *NodalProjectionData;
175 
176 typedef PetscErrorCode (*SgsDDNodalStressEval)(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc);
177 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx);
178 typedef struct {
179   DM                        dm_sgs, dm_dd_inputs, dm_dd_outputs;
180   PetscInt                  num_comp_sgs, num_comp_inputs, num_comp_outputs;
181   OperatorApplyContext      op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx;
182   CeedVector                sgs_nodal_ceed, grad_velo_ceed;
183   SgsDDNodalStressEval      sgs_nodal_eval;
184   SgsDDNodalStressInference sgs_nodal_inference;
185   void                     *sgs_nodal_inference_ctx;
186   PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx);
187 } *SgsDDData;
188 
189 typedef struct {
190   DM                   dm_dd_training;
191   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
192   PetscScalar          filter_widths[16];
193   OperatorApplyContext op_training_data_calc_ctx;
194   NodalProjectionData  filtered_grad_velo_proj;
195   size_t               training_data_array_dims[2];
196   PetscBool            overwrite_training_data;
197 } *SGS_DD_TrainingData;
198 
199 typedef struct {
200   DM                    dm_filter;
201   PetscInt              num_filtered_fields;
202   CeedInt              *num_field_components;
203   PetscInt              field_prim_state, field_velo_prod;
204   OperatorApplyContext  op_rhs_ctx;
205   KSP                   ksp;
206   PetscObjectState      X_loc_state;
207   PetscBool             do_mms_test;
208   CeedContextFieldLabel filter_width_scaling_label;
209 } *DiffFilterData;
210 
211 typedef struct {
212   void    *client;
213   char     rank_id_name[16];
214   PetscInt collocated_database_num_ranks;
215 } *SmartSimData;
216 
217 // PETSc user data
218 struct User_private {
219   MPI_Comm             comm;
220   DM                   dm;
221   DM                   dm_viz;
222   Mat                  interp_viz;
223   Ceed                 ceed;
224   Units                units;
225   Vec                  Q_loc, Q_dot_loc;
226   Physics              phys;
227   AppCtx               app_ctx;
228   CeedVector           q_ceed, q_dot_ceed, g_ceed, x_ceed;
229   CeedOperator         op_ifunction;
230   Mat                  mat_ijacobian;
231   KSP                  mass_ksp;
232   OperatorApplyContext op_rhs_ctx, op_strong_bc_ctx;
233   CeedScalar           time_bc_set;
234   SpanStatsData        spanstats;
235   NodalProjectionData  grad_velo_proj;
236   SgsDDData            sgs_dd_data;
237   DiffFilterData       diff_filter;
238   SmartSimData         smartsim;
239   SGS_DD_TrainingData  sgs_dd_train;
240 };
241 
242 // Units
243 struct Units_private {
244   // fundamental units
245   PetscScalar meter;
246   PetscScalar kilogram;
247   PetscScalar second;
248   PetscScalar Kelvin;
249   // derived units
250   PetscScalar Pascal;
251   PetscScalar J_per_kg_K;
252   PetscScalar m_per_squared_s;
253   PetscScalar W_per_m_K;
254   PetscScalar Joule;
255 };
256 
257 // Boundary conditions
258 struct SimpleBC_private {
259   PetscInt num_inflow, num_outflow, num_freestream, num_slip;
260   PetscInt inflows[16], outflows[16], freestreams[16], slips[16];
261 };
262 
263 // Struct that contains all enums and structs used for the physics of all problems
264 struct Physics_private {
265   PetscBool             implicit;
266   StateVariable         state_var;
267   CeedContextFieldLabel solution_time_label;
268   CeedContextFieldLabel stg_solution_time_label;
269   CeedContextFieldLabel timestep_size_label;
270   CeedContextFieldLabel ics_time_label;
271 };
272 
273 PetscErrorCode BoundaryConditionSetUp(User user, ProblemData problem, AppCtx app_ctx, SimpleBC bc);
274 
275 typedef struct {
276   CeedQFunctionUser    qfunction;
277   const char          *qfunction_loc;
278   CeedQFunctionContext qfunction_context;
279 } ProblemQFunctionSpec;
280 
281 // Problem specific data
282 struct ProblemData_private {
283   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
284   CeedScalar           dm_scale;
285   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
286       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
287   bool          compute_exact_solution_error;
288   PetscBool     set_bc_from_ics, use_strong_bc_ceed, uses_newtonian;
289   size_t        num_bc_defs;
290   BCDefinition *bc_defs;
291   PetscErrorCode (*print_info)(User, ProblemData, AppCtx);
292   PetscErrorCode (*create_mass_operator)(User, CeedOperator *);
293 };
294 
295 extern int FreeContextPetsc(void *);
296 
297 // -----------------------------------------------------------------------------
298 // Set up problems
299 // -----------------------------------------------------------------------------
300 // Set up function for each problem
301 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
302 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
303 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
304 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
305 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
306 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
307 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
308 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
309 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
310 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
311 
312 // Print function for each problem
313 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx);
314 
315 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx);
316 
317 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx);
318 
319 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx);
320 
321 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx);
322 
323 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts);
324 
325 // -----------------------------------------------------------------------------
326 // libCEED functions
327 // -----------------------------------------------------------------------------
328 // Utility function to create local CEED restriction
329 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
330                                          CeedElemRestriction *elem_restr);
331 
332 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
333                                                CeedElemRestriction *restriction);
334 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
335                                                          CeedElemRestriction *restriction);
336 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
337                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
338 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
339                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
340 
341 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
342 PetscErrorCode DMPlexCeedBasisCellToFaceCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face,
343                                                          CeedBasis *basis);
344 PetscErrorCode DMPlexCeedBasisCellToFaceCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt face, PetscInt dm_field,
345                                                CeedBasis *basis);
346 PetscErrorCode DMPlexCreateFaceLabel(DM dm, PetscInt dm_face, PetscInt *num_label_values, char **face_label_name);
347 
348 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc);
349 
350 PetscErrorCode QDataGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
351                         CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
352 PetscErrorCode QDataGetNumComponents(DM dm, CeedInt *q_data_size);
353 PetscErrorCode QDataBoundaryGet(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, CeedElemRestriction elem_restr_x, CeedBasis basis_x,
354                                 CeedVector x_coord, CeedElemRestriction *elem_restr_qd, CeedVector *q_data, CeedInt *q_data_size);
355 PetscErrorCode QDataBoundaryGetNumComponents(DM dm, CeedInt *q_data_size);
356 // -----------------------------------------------------------------------------
357 // Time-stepping functions
358 // -----------------------------------------------------------------------------
359 // RHS (Explicit time-stepper) function setup
360 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data);
361 
362 // Implicit time-stepper function setup
363 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data);
364 
365 // User provided TS Monitor
366 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx);
367 
368 // TS: Create, setup, and solve
369 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, ProblemData problem, Vec *Q, PetscScalar *f_time, TS *ts);
370 
371 // Update Boundary Values when time has changed
372 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t);
373 
374 // -----------------------------------------------------------------------------
375 // Setup DM
376 // -----------------------------------------------------------------------------
377 // Create mesh
378 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData problem, MatType, VecType, DM *dm);
379 
380 // Set up DM
381 PetscErrorCode SetUpDM(DM dm, ProblemData problem, PetscInt degree, PetscInt q_extra, SimpleBC bc, Physics phys);
382 PetscErrorCode DMSetupByOrderBegin_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
383                                        PetscInt num_fields, const PetscInt *field_sizes, DM dm);
384 PetscErrorCode DMSetupByOrderEnd_FEM(PetscBool setup_coords, DM dm);
385 PetscErrorCode DMSetupByOrder_FEM(PetscBool setup_faces, PetscBool setup_coords, PetscInt degree, PetscInt coord_order, PetscInt q_extra,
386                                   PetscInt num_fields, const PetscInt *field_sizes, DM dm);
387 
388 // Refine DM for high-order viz
389 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData problem, SimpleBC bc, Physics phys);
390 
391 // -----------------------------------------------------------------------------
392 // Process command line options
393 // -----------------------------------------------------------------------------
394 // Register problems to be available on the command line
395 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx);
396 
397 // Process general command line options
398 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc);
399 
400 // -----------------------------------------------------------------------------
401 // Miscellaneous utility functions
402 // -----------------------------------------------------------------------------
403 PetscErrorCode GetInverseMultiplicity(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
404                                       PetscBool get_global_multiplicity, CeedElemRestriction *elem_restr_inv_multiplicity,
405                                       CeedVector *inv_multiplicity);
406 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time);
407 
408 PetscErrorCode DMPlexInsertBoundaryValues_FromICs(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM,
409                                                   Vec grad_FVM);
410 
411 // Compare reference solution values with current test run for CI
412 PetscErrorCode RegressionTest(AppCtx app_ctx, Vec Q);
413 
414 // Get error for problems with exact solutions
415 PetscErrorCode PrintError(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time);
416 
417 // Post-processing
418 PetscErrorCode PostProcess(TS ts, CeedData ceed_data, DM dm, ProblemData problem, User user, Vec Q, PetscScalar final_time);
419 
420 // -- Gather initial Q values in case of continuation of simulation
421 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q);
422 
423 // Record boundary values from initial condition
424 PetscErrorCode SetBCsFromICs(DM dm, Vec Q, Vec Q_loc);
425 
426 // Versioning token for binary checkpoints
427 extern const PetscInt32 FLUIDS_FILE_TOKEN;  // for backwards compatibility
428 extern const PetscInt32 FLUIDS_FILE_TOKEN_32;
429 extern const PetscInt32 FLUIDS_FILE_TOKEN_64;
430 
431 // Create appropriate mass qfunction based on number of components N
432 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf);
433 
434 PetscErrorCode NodalProjectionDataDestroy(NodalProjectionData context);
435 
436 PetscErrorCode PhastaDatFileOpen(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], const PetscInt char_array_len, PetscInt dims[2],
437                                  FILE **fp);
438 
439 PetscErrorCode PhastaDatFileGetNRows(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscInt *nrows);
440 
441 PetscErrorCode PhastaDatFileReadToArrayReal(const MPI_Comm comm, const char path[PETSC_MAX_PATH_LEN], PetscReal array[]);
442 
443 // -----------------------------------------------------------------------------
444 // Turbulence Statistics Collection Functions
445 // -----------------------------------------------------------------------------
446 
447 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
448 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
449 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data);
450 
451 // -----------------------------------------------------------------------------
452 // Data-Driven Subgrid Stress (DD-SGS) Modeling Functions
453 // -----------------------------------------------------------------------------
454 
455 PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
456 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data);
457 PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc);
458 PetscErrorCode VelocityGradientProjectionSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, StateVariable state_var_input,
459                                                CeedElemRestriction elem_restr_input, CeedBasis basis_input, NodalProjectionData *pgrad_velo_proj);
460 PetscErrorCode VelocityGradientProjectionApply(NodalProjectionData grad_velo_proj, Vec Q_loc, Vec VelocityGradient);
461 PetscErrorCode GridAnisotropyTensorProjectionSetupApply(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
462                                                         CeedVector *grid_aniso_vector);
463 PetscErrorCode GridAnisotropyTensorCalculateCollocatedVector(Ceed ceed, User user, CeedData ceed_data, CeedElemRestriction *elem_restr_grid_aniso,
464                                                              CeedVector *aniso_colloc_ceed, PetscInt *num_comp_aniso);
465 
466 // -----------------------------------------------------------------------------
467 // Boundary Condition Related Functions
468 // -----------------------------------------------------------------------------
469 
470 // Setup StrongBCs that use QFunctions
471 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc);
472 
473 PetscErrorCode FreestreamBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
474 PetscErrorCode OutflowBCSetup(ProblemData problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference);
475 PetscErrorCode SlipBCSetup(ProblemData problem, DM dm, void *ctx, CeedQFunctionContext newtonian_ig_qfctx);
476 
477 // -----------------------------------------------------------------------------
478 // Differential Filtering Functions
479 // -----------------------------------------------------------------------------
480 
481 PetscErrorCode DifferentialFilterSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
482 PetscErrorCode DifferentialFilterDataDestroy(DiffFilterData diff_filter);
483 PetscErrorCode TSMonitor_DifferentialFilter(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx);
484 PetscErrorCode DifferentialFilterApply(User user, const PetscReal solution_time, const Vec Q, Vec Filtered_Solution);
485 PetscErrorCode DifferentialFilterMmsICSetup(ProblemData problem);
486 
487 // -----------------------------------------------------------------------------
488 // SGS Data-Driven Training via SmartSim
489 // -----------------------------------------------------------------------------
490 PetscErrorCode SmartSimSetup(User user);
491 PetscErrorCode SmartSimDataDestroy(SmartSimData smartsim);
492 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem);
493 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx);
494 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts);
495 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train);
496