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