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