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 // libCEED + PETSc Example: Navier-Stokes 9 // 10 // This example demonstrates a simple usage of libCEED with PETSc to solve a Navier-Stokes problem. 11 // 12 // Build with: 13 // 14 // make [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] navierstokes 15 // 16 // Sample runs: 17 // 18 // ./navierstokes -ceed /cpu/self -problem density_current -degree 1 19 // ./navierstokes -ceed /gpu/cuda -problem advection -degree 1 20 // 21 //TESTARGS(name="blasius") -ceed {ceed_resource} -test -options_file examples/fluids/tests-output/blasius_test.yaml -compare_final_state_atol 2E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius.bin 22 //TESTARGS(name="blasius_STG") -ceed {ceed_resource} -test -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 2E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG.bin 23 //TESTARGS(name="blasius_STG_weakT") -ceed {ceed_resource} -test -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG_weakT.bin -weakT 24 //TESTARGS(name="blasius_STG_strongBC") -ceed {ceed_resource} -test -options_file examples/fluids/tests-output/blasius_stgtest.yaml -compare_final_state_atol 1E-10 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-blasius_STG_strongBC.bin -stg_strong true 25 //TESTARGS(name="channel") -ceed {ceed_resource} -test -options_file examples/fluids/channel.yaml -compare_final_state_atol 2e-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-channel.bin -dm_plex_box_faces 5,5,1 -ts_max_steps 5 26 //TESTARGS(name="channel-primitive") -ceed {ceed_resource} -test -options_file examples/fluids/channel.yaml -compare_final_state_atol 2e-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-channel-prim.bin -dm_plex_box_faces 5,5,1 -ts_max_steps 5 -state_var primitive 27 //TESTARGS(name="dc_explicit") -ceed {ceed_resource} -test -degree 3 -q_extra 2 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -units_kilogram 1e-9 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -mu 75 -ts_dt 1e-3 -units_meter 1e-2 -units_second 1e-2 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-explicit.bin 28 //TESTARGS(name="dc_implicit_stab_none") -ceed {ceed_resource} -test -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -units_kilogram 1e-9 -center 62.5,62.5,187.5 -rc 100. -thetaC -35. -mu 75 -units_meter 1e-2 -units_second 1e-2 -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-dc-implicit-stab-none.bin 29 //TESTARGS(name="adv_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -bc_wall 1,2,3,4,5,6 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-rotation-implicit-stab-supg.bin 30 //TESTARGS(name="adv_translation_implicit_stab_su") -ceed {ceed_resource} -test -problem advection -CtauS .3 -stab su -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -wind_type translation -wind_translation .53,-1.33,-2.65 -bc_inflow 1,2,3,4,5,6 -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv-translation-implicit-stab-su.bin 31 //TESTARGS(name="adv2d_rotation_explicit_strong") -ceed {ceed_resource} -test -problem advection2d -strong_form 1 -degree 3 -dm_plex_box_faces 2,2 -dm_plex_box_lower 0,0 -dm_plex_box_upper 125,125 -bc_wall 1,2,3,4 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ts_dt 1e-3 -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-explicit-strong.bin 32 //TESTARGS(name="adv2d_rotation_implicit_stab_supg") -ceed {ceed_resource} -test -problem advection2d -CtauS .3 -stab supg -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0 -dm_plex_box_upper 125,125 -bc_wall 1,2,3,4 -wall_comps 4 -units_kilogram 1e-9 -rc 100. -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-adv2d-rotation-implicit-stab-supg.bin 33 //TESTARGS(name="euler_implicit") -ceed {ceed_resource} -test -problem euler_vortex -degree 3 -dm_plex_box_faces 1,1,2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_meter 1e-4 -units_second 1e-4 -mean_velocity 1.4,-2.,0 -bc_inflow 4,6 -bc_outflow 3,5 -bc_slip_z 1,2 -vortex_strength 2 -ksp_atol 1e-4 -ksp_rtol 1e-3 -ksp_type bcgs -snes_atol 1e-3 -snes_lag_jacobian 100 -snes_lag_jacobian_persists -snes_mf_operator -ts_dt 1e-3 -implicit -dm_mat_preallocate_skip 0 -ts_type alpha -compare_final_state_atol 5E-4 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-euler-implicit.bin 34 //TESTARGS(name="euler_explicit") -ceed {ceed_resource} -test -problem euler_vortex -degree 3 -q_extra 2 -dm_plex_box_faces 2,2,1 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 125,125,250 -dm_plex_dim 3 -units_meter 1e-4 -units_second 1e-4 -mean_velocity 1.4,-2.,0 -bc_inflow 4,6 -bc_outflow 3,5 -bc_slip_z 1,2 -vortex_strength 2 -ts_dt 1e-7 -ts_rk_type 5bs -ts_rtol 1e-10 -ts_atol 1e-10 -compare_final_state_atol 1E-7 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-euler-explicit.bin 35 //TESTARGS(name="shocktube_explicit_su_yzb") -ceed {ceed_resource} -test -problem shocktube -degree 1 -q_extra 2 -dm_plex_box_faces 50,1,1 -units_meter 1e-2 units_second 1e-2 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1000,20,20 -dm_plex_dim 3 -bc_slip_x 5,6 -bc_slip_y 3,4 -bc_Slip_z 1,2 -yzb -stab su -compare_final_state_atol 1E-11 -compare_final_state_filename examples/fluids/tests-output/fluids-navierstokes-shocktube-explicit-su-yzb.bin 36 37 /// @file 38 /// Navier-Stokes example using PETSc 39 40 const char help[] = "Solve Navier-Stokes using PETSc and libCEED\n"; 41 42 #include "navierstokes.h" 43 44 int main(int argc, char **argv) { 45 // --------------------------------------------------------------------------- 46 // Initialize PETSc 47 // --------------------------------------------------------------------------- 48 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 49 50 // --------------------------------------------------------------------------- 51 // Create structs 52 // --------------------------------------------------------------------------- 53 AppCtx app_ctx; 54 PetscCall(PetscCalloc1(1, &app_ctx)); 55 56 ProblemData *problem = NULL; 57 PetscCall(PetscCalloc1(1, &problem)); 58 59 User user; 60 PetscCall(PetscCalloc1(1, &user)); 61 62 CeedData ceed_data; 63 PetscCall(PetscCalloc1(1, &ceed_data)); 64 65 SimpleBC bc; 66 PetscCall(PetscCalloc1(1, &bc)); 67 68 Physics phys_ctx; 69 PetscCall(PetscCalloc1(1, &phys_ctx)); 70 71 Units units; 72 PetscCall(PetscCalloc1(1, &units)); 73 74 user->app_ctx = app_ctx; 75 user->units = units; 76 user->phys = phys_ctx; 77 problem->bc_from_ics = PETSC_TRUE; 78 79 // --------------------------------------------------------------------------- 80 // Process command line options 81 // --------------------------------------------------------------------------- 82 // -- Register problems to be available on the command line 83 PetscCall(RegisterProblems_NS(app_ctx)); 84 85 // -- Process general command line options 86 MPI_Comm comm = PETSC_COMM_WORLD; 87 user->comm = comm; 88 PetscCall(ProcessCommandLineOptions(comm, app_ctx, bc)); 89 90 // --------------------------------------------------------------------------- 91 // Initialize libCEED 92 // --------------------------------------------------------------------------- 93 // -- Initialize backend 94 Ceed ceed; 95 CeedInit(app_ctx->ceed_resource, &ceed); 96 user->ceed = ceed; 97 98 // -- Check preferred MemType 99 CeedMemType mem_type_backend; 100 CeedGetPreferredMemType(ceed, &mem_type_backend); 101 102 // --------------------------------------------------------------------------- 103 // Set up global mesh 104 // --------------------------------------------------------------------------- 105 // -- Create DM 106 DM dm; 107 VecType vec_type = NULL; 108 MatType mat_type = NULL; 109 switch (mem_type_backend) { 110 case CEED_MEM_HOST: 111 vec_type = VECSTANDARD; 112 break; 113 case CEED_MEM_DEVICE: { 114 const char *resolved; 115 CeedGetResource(ceed, &resolved); 116 if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 117 else if (strstr(resolved, "/gpu/hip")) vec_type = VECKOKKOS; 118 else vec_type = VECSTANDARD; 119 } 120 } 121 if (strstr(vec_type, VECCUDA)) mat_type = MATAIJCUSPARSE; 122 else if (strstr(vec_type, VECKOKKOS)) mat_type = MATAIJKOKKOS; 123 else mat_type = MATAIJ; 124 PetscCall(CreateDM(comm, problem, mat_type, vec_type, &dm)); 125 user->dm = dm; 126 PetscCall(DMSetApplicationContext(dm, user)); 127 128 // --------------------------------------------------------------------------- 129 // Choose the problem from the list of registered problems 130 // --------------------------------------------------------------------------- 131 { 132 PetscErrorCode (*p)(ProblemData *, DM, void *, SimpleBC); 133 PetscCall(PetscFunctionListFind(app_ctx->problems, app_ctx->problem_name, &p)); 134 if (!p) SETERRQ(PETSC_COMM_SELF, 1, "Problem '%s' not found", app_ctx->problem_name); 135 PetscCall((*p)(problem, dm, &user, bc)); 136 } 137 138 // -- Set up DM 139 PetscCall(SetUpDM(dm, problem, app_ctx->degree, bc, phys_ctx)); 140 141 // -- Refine DM for high-order viz 142 if (app_ctx->viz_refine) { 143 PetscCall(VizRefineDM(dm, user, problem, bc, phys_ctx)); 144 } 145 146 // --------------------------------------------------------------------------- 147 // Set up libCEED 148 // --------------------------------------------------------------------------- 149 // -- Set up libCEED objects 150 PetscCall(SetupLibceed(ceed, ceed_data, dm, user, app_ctx, problem, bc)); 151 152 // --------------------------------------------------------------------------- 153 // Set up ICs 154 // --------------------------------------------------------------------------- 155 // -- Set up global state vector Q 156 Vec Q; 157 PetscCall(DMCreateGlobalVector(dm, &Q)); 158 PetscCall(VecZeroEntries(Q)); 159 160 // -- Set up local state vectors Q_loc, Q_dot_loc 161 PetscCall(DMCreateLocalVector(dm, &user->Q_loc)); 162 PetscCall(DMCreateLocalVector(dm, &user->Q_dot_loc)); 163 PetscCall(VecZeroEntries(user->Q_dot_loc)); 164 165 // -- Fix multiplicity for ICs 166 PetscCall(ICs_FixMultiplicity(dm, ceed_data, user, user->Q_loc, Q, 0.0)); 167 168 // --------------------------------------------------------------------------- 169 // Set up lumped mass matrix 170 // --------------------------------------------------------------------------- 171 // -- Set up global mass vector 172 PetscCall(VecDuplicate(Q, &user->M)); 173 174 // -- Compute lumped mass matrix 175 PetscCall(ComputeLumpedMassMatrix(ceed, dm, ceed_data, user->M)); 176 177 // --------------------------------------------------------------------------- 178 // Record boundary values from initial condition 179 // --------------------------------------------------------------------------- 180 // -- This overrides DMPlexInsertBoundaryValues(). 181 // We use this for the main simulation DM because the reference DMPlexInsertBoundaryValues() is very slow on the GPU due to extra device-to-host 182 // communication. If we disable this, we should still get the same results due to the problem->bc function, but with potentially much slower 183 // execution. 184 if (problem->bc_from_ics) { 185 PetscCall(SetBCsFromICs_NS(dm, Q, user->Q_loc)); 186 } 187 188 // --------------------------------------------------------------------------- 189 // Create output directory 190 // --------------------------------------------------------------------------- 191 PetscMPIInt rank; 192 MPI_Comm_rank(comm, &rank); 193 if (!rank) { 194 PetscCall(PetscMkdir(app_ctx->output_dir)); 195 } 196 197 // --------------------------------------------------------------------------- 198 // Gather initial Q values in case of continuation of simulation 199 // --------------------------------------------------------------------------- 200 // -- Set up initial values from binary file 201 if (app_ctx->cont_steps) { 202 PetscCall(SetupICsFromBinary(comm, app_ctx, Q)); 203 } 204 205 // --------------------------------------------------------------------------- 206 // Print problem summary 207 // --------------------------------------------------------------------------- 208 if (!app_ctx->test_mode) { 209 // Header and rank 210 char host_name[PETSC_MAX_PATH_LEN]; 211 int comm_size; 212 PetscCall(PetscGetHostName(host_name, sizeof host_name)); 213 PetscCall(MPI_Comm_size(comm, &comm_size)); 214 PetscCall(PetscPrintf(comm, 215 "\n-- Navier-Stokes solver - libCEED + PETSc --\n" 216 " MPI:\n" 217 " Host Name : %s\n" 218 " Total ranks : %d\n", 219 host_name, comm_size)); 220 221 // Problem specific info 222 PetscCall(problem->print_info(problem, app_ctx)); 223 224 // libCEED 225 const char *used_resource; 226 CeedGetResource(ceed, &used_resource); 227 PetscCall(PetscPrintf(comm, 228 " libCEED:\n" 229 " libCEED Backend : %s\n" 230 " libCEED Backend MemType : %s\n", 231 used_resource, CeedMemTypes[mem_type_backend])); 232 // PETSc 233 char box_faces_str[PETSC_MAX_PATH_LEN] = "3,3,3"; 234 if (problem->dim == 2) box_faces_str[3] = '\0'; 235 PetscCall(PetscOptionsGetString(NULL, NULL, "-dm_plex_box_faces", box_faces_str, sizeof(box_faces_str), NULL)); 236 MatType mat_type; 237 VecType vec_type; 238 PetscCall(DMGetMatType(dm, &mat_type)); 239 PetscCall(DMGetVecType(dm, &vec_type)); 240 PetscCall(PetscPrintf(comm, 241 " PETSc:\n" 242 " Box Faces : %s\n" 243 " DM MatType : %s\n" 244 " DM VecType : %s\n" 245 " Time Stepping Scheme : %s\n", 246 box_faces_str, mat_type, vec_type, phys_ctx->implicit ? "implicit" : "explicit")); 247 // Mesh 248 const PetscInt num_comp_q = 5; 249 CeedInt glob_dofs, owned_dofs; 250 PetscInt glob_nodes, local_nodes; 251 const CeedInt num_P = app_ctx->degree + 1, num_Q = num_P + app_ctx->q_extra; 252 // -- Get global size 253 PetscCall(VecGetSize(Q, &glob_dofs)); 254 PetscCall(VecGetLocalSize(Q, &owned_dofs)); 255 glob_nodes = glob_dofs / num_comp_q; 256 // -- Get local size 257 PetscCall(VecGetSize(user->Q_loc, &local_nodes)); 258 local_nodes /= num_comp_q; 259 PetscCall(PetscPrintf(comm, 260 " Mesh:\n" 261 " Number of 1D Basis Nodes (P) : %" CeedInt_FMT "\n" 262 " Number of 1D Quadrature Points (Q) : %" CeedInt_FMT "\n" 263 " Global DoFs : %" PetscInt_FMT "\n" 264 " Owned DoFs : %" PetscInt_FMT "\n" 265 " DoFs per node : %" PetscInt_FMT "\n" 266 " Global nodes (DoFs / %" PetscInt_FMT ") : %" PetscInt_FMT "\n" 267 " Local nodes : %" PetscInt_FMT "\n", 268 num_P, num_Q, glob_dofs, owned_dofs, num_comp_q, num_comp_q, glob_nodes, local_nodes)); 269 } 270 // -- Zero Q_loc 271 PetscCall(VecZeroEntries(user->Q_loc)); 272 273 // --------------------------------------------------------------------------- 274 // TS: Create, setup, and solve 275 // --------------------------------------------------------------------------- 276 TS ts; 277 PetscScalar final_time; 278 PetscCall(TSSolve_NS(dm, user, app_ctx, phys_ctx, &Q, &final_time, &ts)); 279 280 // --------------------------------------------------------------------------- 281 // Post-processing 282 // --------------------------------------------------------------------------- 283 PetscCall(PostProcess_NS(ts, ceed_data, dm, problem, user, Q, final_time)); 284 285 // --------------------------------------------------------------------------- 286 // Destroy libCEED objects 287 // --------------------------------------------------------------------------- 288 // -- Vectors 289 CeedVectorDestroy(&ceed_data->x_coord); 290 CeedVectorDestroy(&ceed_data->q_data); 291 CeedVectorDestroy(&user->q_ceed); 292 CeedVectorDestroy(&user->q_dot_ceed); 293 CeedVectorDestroy(&user->g_ceed); 294 CeedVectorDestroy(&user->coo_values_amat); 295 CeedVectorDestroy(&user->coo_values_pmat); 296 297 // -- QFunctions 298 CeedQFunctionDestroy(&ceed_data->qf_setup_vol); 299 CeedQFunctionDestroy(&ceed_data->qf_ics); 300 CeedQFunctionDestroy(&ceed_data->qf_rhs_vol); 301 CeedQFunctionDestroy(&ceed_data->qf_ifunction_vol); 302 CeedQFunctionDestroy(&ceed_data->qf_setup_sur); 303 CeedQFunctionDestroy(&ceed_data->qf_apply_inflow); 304 CeedQFunctionDestroy(&ceed_data->qf_apply_inflow_jacobian); 305 CeedQFunctionDestroy(&ceed_data->qf_apply_outflow); 306 CeedQFunctionDestroy(&ceed_data->qf_apply_outflow_jacobian); 307 308 // -- Bases 309 CeedBasisDestroy(&ceed_data->basis_q); 310 CeedBasisDestroy(&ceed_data->basis_x); 311 CeedBasisDestroy(&ceed_data->basis_xc); 312 CeedBasisDestroy(&ceed_data->basis_q_sur); 313 CeedBasisDestroy(&ceed_data->basis_x_sur); 314 315 // -- Restrictions 316 CeedElemRestrictionDestroy(&ceed_data->elem_restr_q); 317 CeedElemRestrictionDestroy(&ceed_data->elem_restr_x); 318 CeedElemRestrictionDestroy(&ceed_data->elem_restr_qd_i); 319 320 // -- Operators 321 CeedOperatorDestroy(&ceed_data->op_setup_vol); 322 CeedOperatorDestroy(&ceed_data->op_ics); 323 CeedOperatorDestroy(&user->op_rhs_vol); 324 CeedOperatorDestroy(&user->op_ifunction_vol); 325 CeedOperatorDestroy(&user->op_rhs); 326 CeedOperatorDestroy(&user->op_ifunction); 327 CeedOperatorDestroy(&user->op_ijacobian); 328 329 // -- Ceed 330 CeedDestroy(&ceed); 331 332 // --------------------------------------------------------------------------- 333 // Clean up PETSc 334 // --------------------------------------------------------------------------- 335 // -- Vectors 336 PetscCall(VecDestroy(&Q)); 337 PetscCall(VecDestroy(&user->M)); 338 PetscCall(VecDestroy(&user->Q_loc)); 339 PetscCall(VecDestroy(&user->Q_dot_loc)); 340 341 // -- Matrices 342 PetscCall(MatDestroy(&user->interp_viz)); 343 344 // -- DM 345 PetscCall(DMDestroy(&dm)); 346 PetscCall(DMDestroy(&user->dm_viz)); 347 348 // -- TS 349 PetscCall(TSDestroy(&ts)); 350 351 // -- Function list 352 PetscCall(PetscFunctionListDestroy(&app_ctx->problems)); 353 354 PetscCall(PetscFree(app_ctx->amat_type)); 355 356 // -- Structs 357 PetscCall(PetscFree(units)); 358 PetscCall(PetscFree(user)); 359 PetscCall(PetscFree(problem)); 360 PetscCall(PetscFree(bc)); 361 PetscCall(PetscFree(phys_ctx)); 362 PetscCall(PetscFree(app_ctx)); 363 PetscCall(PetscFree(ceed_data)); 364 365 return PetscFinalize(); 366 } 367