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