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