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 <petscdm.h> 13 #include <petscdmplex.h> 14 #include <petscsys.h> 15 #include <petscts.h> 16 #include <stdbool.h> 17 18 #include "./include/matops.h" 19 #include "qfunctions/newtonian_types.h" 20 #include "qfunctions/stabilization_types.h" 21 22 // ----------------------------------------------------------------------------- 23 // PETSc Version 24 // ----------------------------------------------------------------------------- 25 #if PETSC_VERSION_LT(3, 17, 0) 26 #error "PETSc v3.17 or later is required" 27 #endif 28 29 // ----------------------------------------------------------------------------- 30 // Enums 31 // ----------------------------------------------------------------------------- 32 // Translate PetscMemType to CeedMemType 33 static inline CeedMemType MemTypeP2C(PetscMemType mem_type) { return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST; } 34 35 // Advection - Wind Options 36 typedef enum { 37 WIND_ROTATION = 0, 38 WIND_TRANSLATION = 1, 39 } WindType; 40 static const char *const WindTypes[] = {"rotation", "translation", "WindType", "WIND_", NULL}; 41 42 // Advection - Bubble Types 43 typedef enum { 44 BUBBLE_SPHERE = 0, // dim=3 45 BUBBLE_CYLINDER = 1, // dim=2 46 } BubbleType; 47 static const char *const BubbleTypes[] = {"sphere", "cylinder", "BubbleType", "BUBBLE_", NULL}; 48 49 // Advection - Bubble Continuity Types 50 typedef enum { 51 BUBBLE_CONTINUITY_SMOOTH = 0, // Original continuous, smooth shape 52 BUBBLE_CONTINUITY_BACK_SHARP = 1, // Discontinuous, sharp back half shape 53 BUBBLE_CONTINUITY_THICK = 2, // Define a finite thickness 54 } BubbleContinuityType; 55 static const char *const BubbleContinuityTypes[] = {"smooth", "back_sharp", "thick", "BubbleContinuityType", "BUBBLE_CONTINUITY_", NULL}; 56 57 // Euler - test cases 58 typedef enum { 59 EULER_TEST_ISENTROPIC_VORTEX = 0, 60 EULER_TEST_1 = 1, 61 EULER_TEST_2 = 2, 62 EULER_TEST_3 = 3, 63 EULER_TEST_4 = 4, 64 EULER_TEST_5 = 5, 65 } EulerTestType; 66 static const char *const EulerTestTypes[] = {"isentropic_vortex", "test_1", "test_2", "test_3", "test_4", "test_5", 67 "EulerTestType", "EULER_TEST_", NULL}; 68 69 // Stabilization methods 70 static const char *const StabilizationTypes[] = {"none", "SU", "SUPG", "StabilizationType", "STAB_", NULL}; 71 72 // Euler - test cases 73 typedef enum { 74 TESTTYPE_NONE = 0, 75 TESTTYPE_SOLVER = 1, 76 TESTTYPE_TURB_SPANSTATS = 2, 77 } TestType; 78 static const char *const TestTypes[] = {"none", "solver", "turb_spanstats", "TestType", "TESTTYPE_", NULL}; 79 80 // ----------------------------------------------------------------------------- 81 // Structs 82 // ----------------------------------------------------------------------------- 83 // Structs declarations 84 typedef struct AppCtx_private *AppCtx; 85 typedef struct CeedData_private *CeedData; 86 typedef struct User_private *User; 87 typedef struct Units_private *Units; 88 typedef struct SimpleBC_private *SimpleBC; 89 typedef struct Physics_private *Physics; 90 91 // Application context from user command line options 92 struct AppCtx_private { 93 // libCEED arguments 94 char ceed_resource[PETSC_MAX_PATH_LEN]; // libCEED backend 95 PetscInt degree; 96 PetscInt q_extra; 97 // Solver arguments 98 MatType amat_type; 99 PetscBool pmat_pbdiagonal; 100 // Post-processing arguments 101 PetscInt checkpoint_interval; 102 PetscInt viz_refine; 103 PetscInt cont_steps; 104 PetscReal cont_time; 105 char cont_file[PETSC_MAX_PATH_LEN]; 106 char cont_time_file[PETSC_MAX_PATH_LEN]; 107 char output_dir[PETSC_MAX_PATH_LEN]; 108 PetscBool add_stepnum2bin; 109 PetscBool checkpoint_vtk; 110 // Problem type arguments 111 PetscFunctionList problems; 112 char problem_name[PETSC_MAX_PATH_LEN]; 113 // Test mode arguments 114 TestType test_type; 115 PetscScalar test_tol; 116 char test_file_path[PETSC_MAX_PATH_LEN]; 117 // Turbulent spanwise statistics 118 PetscBool turb_spanstats_enable; 119 PetscInt turb_spanstats_collect_interval; 120 PetscInt turb_spanstats_viewer_interval; 121 PetscViewer turb_spanstats_viewer; 122 PetscViewerFormat turb_spanstats_viewer_format; 123 // Wall forces 124 struct { 125 PetscInt num_wall; 126 PetscInt *walls; 127 PetscViewer viewer; 128 PetscViewerFormat viewer_format; 129 PetscBool header_written; 130 } wall_forces; 131 }; 132 133 // libCEED data struct 134 struct CeedData_private { 135 CeedVector x_coord, q_data; 136 CeedBasis basis_x, basis_xc, basis_q, basis_x_sur, basis_q_sur, basis_xc_sur; 137 CeedElemRestriction elem_restr_x, elem_restr_q, elem_restr_qd_i; 138 CeedOperator op_setup_vol, op_ics; 139 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, 140 qf_apply_outflow_jacobian, qf_apply_freestream, qf_apply_freestream_jacobian; 141 }; 142 143 typedef struct { 144 DM dm; 145 PetscSF sf; // For communicating child data to parents 146 CeedOperator op_stats_collect, op_stats_proj; 147 PetscInt num_comp_stats; 148 CeedVector child_stats, parent_stats; // collocated statistics data 149 CeedVector rhs_ceed; 150 Vec M_inv; // Lumped Mass matrix inverse 151 KSP ksp; // For the L^2 projection solve 152 CeedScalar span_width; // spanwise width of the child domain 153 PetscBool do_mms_test; 154 MatopApplyContext mms_error_ctx; 155 CeedContextFieldLabel solution_time_label, previous_time_label; 156 } Span_Stats; 157 158 // PETSc user data 159 struct User_private { 160 MPI_Comm comm; 161 DM dm; 162 DM dm_viz; 163 Mat interp_viz; 164 Ceed ceed; 165 Units units; 166 Vec M, Q_loc, Q_dot_loc; 167 Physics phys; 168 AppCtx app_ctx; 169 CeedVector q_ceed, q_dot_ceed, g_ceed, coo_values_amat, coo_values_pmat, x_ceed; 170 CeedOperator op_rhs_vol, op_rhs, op_ifunction_vol, op_ifunction, op_ijacobian, op_dirichlet; 171 bool matrices_set_up; 172 CeedScalar time, dt, time_bc_set; 173 Span_Stats spanstats; 174 }; 175 176 // Units 177 struct Units_private { 178 // fundamental units 179 PetscScalar meter; 180 PetscScalar kilogram; 181 PetscScalar second; 182 PetscScalar Kelvin; 183 // derived units 184 PetscScalar Pascal; 185 PetscScalar J_per_kg_K; 186 PetscScalar m_per_squared_s; 187 PetscScalar W_per_m_K; 188 PetscScalar Joule; 189 }; 190 191 // Boundary conditions 192 struct SimpleBC_private { 193 PetscInt num_wall, // Number of faces with wall BCs 194 wall_comps[5], // An array of constrained component numbers 195 num_comps, 196 num_slip[3], // Number of faces with slip BCs 197 num_inflow, num_outflow, num_freestream; 198 PetscInt walls[16], slips[3][16], inflows[16], outflows[16], freestreams[16]; 199 PetscBool user_bc; 200 }; 201 202 // Struct that contains all enums and structs used for the physics of all problems 203 struct Physics_private { 204 WindType wind_type; 205 BubbleType bubble_type; 206 BubbleContinuityType bubble_continuity_type; 207 EulerTestType euler_test; 208 StabilizationType stab; 209 PetscBool implicit; 210 StateVariable state_var; 211 PetscBool has_curr_time; 212 PetscBool has_neumann; 213 CeedContextFieldLabel solution_time_label; 214 CeedContextFieldLabel stg_solution_time_label; 215 CeedContextFieldLabel timestep_size_label; 216 CeedContextFieldLabel ics_time_label; 217 CeedContextFieldLabel ijacobian_time_shift_label; 218 }; 219 220 typedef struct { 221 CeedQFunctionUser qfunction; 222 const char *qfunction_loc; 223 CeedQFunctionContext qfunction_context; 224 } ProblemQFunctionSpec; 225 226 // Problem specific data 227 typedef struct ProblemData_private ProblemData; 228 struct ProblemData_private { 229 CeedInt dim, q_data_size_vol, q_data_size_sur, jac_data_size_sur; 230 CeedScalar dm_scale; 231 ProblemQFunctionSpec setup_vol, setup_sur, ics, apply_vol_rhs, apply_vol_ifunction, apply_vol_ijacobian, apply_inflow, apply_outflow, 232 apply_freestream, apply_inflow_jacobian, apply_outflow_jacobian, apply_freestream_jacobian; 233 bool non_zero_time; 234 PetscErrorCode (*bc)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 235 void *bc_ctx; 236 PetscBool bc_from_ics, use_dirichlet_ceed; 237 PetscErrorCode (*print_info)(ProblemData *, AppCtx); 238 }; 239 240 extern int FreeContextPetsc(void *); 241 242 // ----------------------------------------------------------------------------- 243 // Set up problems 244 // ----------------------------------------------------------------------------- 245 // Set up function for each problem 246 extern PetscErrorCode NS_NEWTONIAN_WAVE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 247 extern PetscErrorCode NS_CHANNEL(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 248 extern PetscErrorCode NS_BLASIUS(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 249 extern PetscErrorCode NS_NEWTONIAN_IG(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 250 extern PetscErrorCode NS_DENSITY_CURRENT(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 251 extern PetscErrorCode NS_EULER_VORTEX(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 252 extern PetscErrorCode NS_SHOCKTUBE(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 253 extern PetscErrorCode NS_ADVECTION(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 254 extern PetscErrorCode NS_ADVECTION2D(ProblemData *problem, DM dm, void *ctx, SimpleBC bc); 255 256 // Print function for each problem 257 extern PetscErrorCode PRINT_NEWTONIAN(ProblemData *problem, AppCtx app_ctx); 258 259 extern PetscErrorCode PRINT_EULER_VORTEX(ProblemData *problem, AppCtx app_ctx); 260 261 extern PetscErrorCode PRINT_SHOCKTUBE(ProblemData *problem, AppCtx app_ctx); 262 263 extern PetscErrorCode PRINT_ADVECTION(ProblemData *problem, AppCtx app_ctx); 264 265 extern PetscErrorCode PRINT_ADVECTION2D(ProblemData *problem, AppCtx app_ctx); 266 267 // ----------------------------------------------------------------------------- 268 // libCEED functions 269 // ----------------------------------------------------------------------------- 270 // Utility function - essential BC dofs are encoded in closure indices as -(i+1). 271 PetscInt Involute(PetscInt i); 272 273 // Utility function to create local CEED restriction 274 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr); 275 276 // Utility function to get Ceed Restriction for each domain 277 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size, 278 CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i); 279 280 // Utility function to create CEED Composite Operator for the entire domain 281 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol, 282 CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur, 283 CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian); 284 285 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc); 286 287 // ----------------------------------------------------------------------------- 288 // Time-stepping functions 289 // ----------------------------------------------------------------------------- 290 // Compute mass matrix for explicit scheme 291 PetscErrorCode ComputeLumpedMassMatrix(Ceed ceed, DM dm, CeedData ceed_data, Vec M); 292 293 // RHS (Explicit time-stepper) function setup 294 PetscErrorCode RHS_NS(TS ts, PetscReal t, Vec Q, Vec G, void *user_data); 295 296 // Implicit time-stepper function setup 297 PetscErrorCode IFunction_NS(TS ts, PetscReal t, Vec Q, Vec Q_dot, Vec G, void *user_data); 298 299 // User provided TS Monitor 300 PetscErrorCode TSMonitor_NS(TS ts, PetscInt step_no, PetscReal time, Vec Q, void *ctx); 301 302 // TS: Create, setup, and solve 303 PetscErrorCode TSSolve_NS(DM dm, User user, AppCtx app_ctx, Physics phys, Vec *Q, PetscScalar *f_time, TS *ts); 304 305 // Update Boundary Values when time has changed 306 PetscErrorCode UpdateBoundaryValues(User user, Vec Q_loc, PetscReal t); 307 308 // ----------------------------------------------------------------------------- 309 // Setup DM 310 // ----------------------------------------------------------------------------- 311 // Create mesh 312 PetscErrorCode CreateDM(MPI_Comm comm, ProblemData *problem, MatType, VecType, DM *dm); 313 314 // Set up DM 315 PetscErrorCode SetUpDM(DM dm, ProblemData *problem, PetscInt degree, SimpleBC bc, Physics phys); 316 317 // Refine DM for high-order viz 318 PetscErrorCode VizRefineDM(DM dm, User user, ProblemData *problem, SimpleBC bc, Physics phys); 319 320 // ----------------------------------------------------------------------------- 321 // Process command line options 322 // ----------------------------------------------------------------------------- 323 // Register problems to be available on the command line 324 PetscErrorCode RegisterProblems_NS(AppCtx app_ctx); 325 326 // Process general command line options 327 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc); 328 329 // ----------------------------------------------------------------------------- 330 // Miscellaneous utility functions 331 // ----------------------------------------------------------------------------- 332 PetscErrorCode ICs_FixMultiplicity(DM dm, CeedData ceed_data, User user, Vec Q_loc, Vec Q, CeedScalar time); 333 334 PetscErrorCode DMPlexInsertBoundaryValues_NS(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, Vec cell_geom_FVM, 335 Vec grad_FVM); 336 337 // Compare reference solution values with current test run for CI 338 PetscErrorCode RegressionTests_NS(AppCtx app_ctx, Vec Q); 339 340 // Get error for problems with exact solutions 341 PetscErrorCode GetError_NS(CeedData ceed_data, DM dm, User user, Vec Q, PetscScalar final_time); 342 343 // Post-processing 344 PetscErrorCode PostProcess_NS(TS ts, CeedData ceed_data, DM dm, ProblemData *problem, User user, Vec Q, PetscScalar final_time); 345 346 // -- Gather initial Q values in case of continuation of simulation 347 PetscErrorCode SetupICsFromBinary(MPI_Comm comm, AppCtx app_ctx, Vec Q); 348 349 // Record boundary values from initial condition 350 PetscErrorCode SetBCsFromICs_NS(DM dm, Vec Q, Vec Q_loc); 351 352 // Versioning token for binary checkpoints 353 extern const PetscInt FLUIDS_FILE_TOKEN; 354 355 // Create appropriate mass qfunction based on number of components N 356 PetscErrorCode CreateMassQFunction(Ceed ceed, CeedInt N, CeedInt q_data_size, CeedQFunction *qf); 357 358 PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc); 359 360 PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem); 361 362 PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx); 363 364 PetscErrorCode DestroyStats(User user, CeedData ceed_data); 365 366 // ----------------------------------------------------------------------------- 367 // Boundary Condition Related Functions 368 // ----------------------------------------------------------------------------- 369 370 // Setup StrongBCs that use QFunctions 371 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur, 372 CeedInt q_data_size_sur); 373 374 PetscErrorCode FreestreamBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 375 PetscErrorCode OutflowBCSetup(ProblemData *problem, DM dm, void *ctx, NewtonianIdealGasContext newtonian_ig_ctx, const StatePrimitive *reference); 376 377 #endif // libceed_fluids_examples_navier_stokes_h 378