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