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