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