xref: /libCEED/examples/fluids/navierstokes.h (revision fe96005463bdbb79b892d21a5c89e2b475ecf62b)
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 <mat-ceed.h>
11 #include <petsc-ceed-utils.h>
12 #include <petscts.h>
13 #include <stdbool.h>
14 
15 #include "./include/petsc_ops.h"
16 #include "qfunctions/newtonian_types.h"
17 
18 #if PETSC_VERSION_LT(3, 21, 0)
19 #error "PETSc v3.21 or later is required"
20 #endif
21 
22 // -----------------------------------------------------------------------------
23 // Enums
24 // -----------------------------------------------------------------------------
25 
26 // Euler - test cases
27 typedef enum {
28   EULER_TEST_ISENTROPIC_VORTEX = 0,
29   EULER_TEST_1                 = 1,
30   EULER_TEST_2                 = 2,
31   EULER_TEST_3                 = 3,
32   EULER_TEST_4                 = 4,
33   EULER_TEST_5                 = 5,
34 } EulerTestType;
35 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1",      "test_2", "test_3", "test_4", "test_5",
36                                              "EulerTestType",     "EULER_TEST_", NULL};
37 
38 // Advection - Wind types
39 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL};
40 
41 // Advection - Initial Condition Types
42 static const char *const AdvectionICTypes[] = {"sphere", "cylinder", "cosine_hill", "skew", "AdvectionICType", "ADVECTIONIC_", NULL};
43 
44 // Advection - Bubble Continuity Types
45 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "cosine", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL};
46 
47 // Stabilization methods
48 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL};
49 
50 // Stabilization tau constants
51 static const char *const StabilizationTauTypes[] = {"Ctau", "AdvDiff_Shakib", "AdvDiff_Shakib_P", "StabilizationTauType", "STAB_TAU_", NULL};
52 
53 // Test mode type
54 typedef enum {
55   TESTTYPE_NONE           = 0,
56   TESTTYPE_SOLVER         = 1,
57   TESTTYPE_TURB_SPANSTATS = 2,
58   TESTTYPE_DIFF_FILTER    = 3,
59 } TestType;
60 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "diff_filter", "TestType", "TESTTYPE_", NULL};
61 
62 // Subgrid-Stress mode type
63 typedef enum {
64   SGS_MODEL_NONE        = 0,
65   SGS_MODEL_DATA_DRIVEN = 1,
66 } SGSModelType;
67 static const char *const SGSModelTypes[] = {"none", "data_driven", "SGSModelType", "SGS_MODEL_", NULL};
68 
69 // Mesh transformation type
70 typedef enum {
71   MESH_TRANSFORM_NONE      = 0,
72   MESH_TRANSFORM_PLATEMESH = 1,
73 } MeshTransformType;
74 static const char *const MeshTransformTypes[] = {"none", "platemesh", "MeshTransformType", "MESH_TRANSFORM_", NULL};
75 
76 static const char *const DifferentialFilterDampingFunctions[] = {
77     "none", "van_driest", "mms", "DifferentialFilterDampingFunction", "DIFF_FILTER_DAMP_", NULL};
78 
79 // -----------------------------------------------------------------------------
80 // Log Events
81 // -----------------------------------------------------------------------------
82 extern PetscLogEvent FLUIDS_CeedOperatorApply;
83 extern PetscLogEvent FLUIDS_SmartRedis_Init;
84 extern PetscLogEvent FLUIDS_SmartRedis_Meta;
85 extern PetscLogEvent FLUIDS_SmartRedis_Train;
86 extern PetscLogEvent FLUIDS_TrainDataCompute;
87 extern PetscLogEvent FLUIDS_DifferentialFilter;
88 extern PetscLogEvent FLUIDS_VelocityGradientProjection;
89 PetscErrorCode       RegisterLogEvents();
90 
91 // -----------------------------------------------------------------------------
92 // Structs
93 // -----------------------------------------------------------------------------
94 // Structs declarations
95 typedef struct AppCtx_private   *AppCtx;
96 typedef struct CeedData_private *CeedData;
97 typedef struct User_private     *User;
98 typedef struct Units_private    *Units;
99 typedef struct SimpleBC_private *SimpleBC;
100 typedef struct Physics_private  *Physics;
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_wall,  // Number of faces with wall BCs
261       wall_comps[5],  // An array of constrained component numbers
262       num_comps,
263       num_symmetry[3],  // Number of faces with symmetry BCs
264       num_inflow, num_outflow, num_freestream, num_slip;
265   PetscInt walls[16], symmetries[3][16], inflows[16], outflows[16], freestreams[16], slips[16];
266 };
267 
268 // Struct that contains all enums and structs used for the physics of all problems
269 struct Physics_private {
270   PetscBool             implicit;
271   StateVariable         state_var;
272   CeedContextFieldLabel solution_time_label;
273   CeedContextFieldLabel stg_solution_time_label;
274   CeedContextFieldLabel timestep_size_label;
275   CeedContextFieldLabel ics_time_label;
276 };
277 
278 typedef struct {
279   CeedQFunctionUser    qfunction;
280   const char          *qfunction_loc;
281   CeedQFunctionContext qfunction_context;
282 } ProblemQFunctionSpec;
283 
284 // Problem specific data
285 typedef struct ProblemData_private *ProblemData;
286 struct ProblemData_private {
287   CeedInt              dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur;
288   CeedScalar           dm_scale;
289   ProblemQFunctionSpec ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, apply_freestream, apply_slip,
290       apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian, apply_slip_jacobian;
291   bool      non_zero_time;
292   PetscBool bc_from_ics, use_strong_bc_ceed, uses_newtonian;
293   PetscErrorCode (*print_info)(User, ProblemData, AppCtx);
294   PetscErrorCode (*create_mass_operator)(User, CeedOperator *);
295 };
296 
297 extern int FreeContextPetsc(void *);
298 
299 // -----------------------------------------------------------------------------
300 // Set up problems
301 // -----------------------------------------------------------------------------
302 // Set up function for each problem
303 extern PetscErrorCode NS_TAYLOR_GREEN(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
304 extern PetscErrorCode NS_GAUSSIAN_WAVE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
305 extern PetscErrorCode NS_CHANNEL(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
306 extern PetscErrorCode NS_BLASIUS(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
307 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
308 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
309 extern PetscErrorCode NS_EULER_VORTEX(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
310 extern PetscErrorCode NS_SHOCKTUBE(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
311 extern PetscErrorCode NS_ADVECTION(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
312 extern PetscErrorCode NS_ADVECTION2D(ProblemData problem, DM dm, void *ctx, SimpleBC bc);
313 
314 // Print function for each problem
315 extern PetscErrorCode PRINT_NEWTONIAN(User user, ProblemData problem, AppCtx app_ctx);
316 
317 extern PetscErrorCode PRINT_EULER_VORTEX(User user, ProblemData problem, AppCtx app_ctx);
318 
319 extern PetscErrorCode PRINT_SHOCKTUBE(User user, ProblemData problem, AppCtx app_ctx);
320 
321 extern PetscErrorCode PRINT_ADVECTION(User user, ProblemData problem, AppCtx app_ctx);
322 
323 extern PetscErrorCode PRINT_ADVECTION2D(User user, ProblemData problem, AppCtx app_ctx);
324 
325 PetscErrorCode PrintRunInfo(User user, Physics phys_ctx, ProblemData problem, TS ts);
326 
327 // -----------------------------------------------------------------------------
328 // libCEED functions
329 // -----------------------------------------------------------------------------
330 // Utility function to create local CEED restriction
331 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
332                                          CeedElemRestriction *elem_restr);
333 
334 PetscErrorCode DMPlexCeedElemRestrictionCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt dm_field,
335                                                CeedElemRestriction *restriction);
336 PetscErrorCode DMPlexCeedElemRestrictionCoordinateCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
337                                                          CeedElemRestriction *restriction);
338 PetscErrorCode DMPlexCeedElemRestrictionQDataCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
339                                                     PetscInt q_data_size, CeedElemRestriction *restriction);
340 PetscErrorCode DMPlexCeedElemRestrictionCollocatedCreate(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height,
341                                                          PetscInt q_data_size, CeedElemRestriction *restriction);
342 
343 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis);
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