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