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