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