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