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