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