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