13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 377841947SLeila Ghaffari // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 577841947SLeila Ghaffari // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 777841947SLeila Ghaffari 877841947SLeila Ghaffari /// @file 977841947SLeila Ghaffari /// Command line option processing for Navier-Stokes example using PETSc 1077841947SLeila Ghaffari 110992c3d1SJed Brown #include <petscdevice.h> 1249aac155SJeremy L Thompson #include <petscsys.h> 130992c3d1SJed Brown 1477841947SLeila Ghaffari #include "../navierstokes.h" 1577841947SLeila Ghaffari 1677841947SLeila Ghaffari // Register problems to be available on the command line 1777841947SLeila Ghaffari PetscErrorCode RegisterProblems_NS(AppCtx app_ctx) { 1877841947SLeila Ghaffari app_ctx->problems = NULL; 1977841947SLeila Ghaffari 20*f17d818dSJames Wright PetscFunctionBeginUser; 212b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "density_current", NS_DENSITY_CURRENT)); 222b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "euler_vortex", NS_EULER_VORTEX)); 232b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "shocktube", NS_SHOCKTUBE)); 242b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection", NS_ADVECTION)); 252b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection2d", NS_ADVECTION2D)); 262b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "blasius", NS_BLASIUS)); 272b730f8bSJeremy L Thompson PetscCall(PetscFunctionListAdd(&app_ctx->problems, "channel", NS_CHANNEL)); 28530ad8c4SKenneth E. Jansen PetscCall(PetscFunctionListAdd(&app_ctx->problems, "gaussian_wave", NS_GAUSSIAN_WAVE)); 29d310b3d3SAdeleke O. Bankole PetscCall(PetscFunctionListAdd(&app_ctx->problems, "newtonian", NS_NEWTONIAN_IG)); 3014712a6bSJames Wright PetscCall(PetscFunctionListAdd(&app_ctx->problems, "taylor_green", NS_TAYLOR_GREEN)); 31ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3277841947SLeila Ghaffari } 3377841947SLeila Ghaffari 3477841947SLeila Ghaffari // Process general command line options 352b730f8bSJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc) { 3677841947SLeila Ghaffari PetscBool ceed_flag = PETSC_FALSE; 3777841947SLeila Ghaffari PetscBool problem_flag = PETSC_FALSE; 3869293791SJames Wright PetscBool option_set = PETSC_FALSE; 392b730f8bSJeremy L Thompson 4077841947SLeila Ghaffari PetscFunctionBeginUser; 412b730f8bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", NULL); 4277841947SLeila Ghaffari 432b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource, 442b730f8bSJeremy L Thompson sizeof(app_ctx->ceed_resource), &ceed_flag)); 4577841947SLeila Ghaffari 468fb33541SJames Wright app_ctx->test_type = TESTTYPE_NONE; 478fb33541SJames Wright PetscCall(PetscOptionsEnum("-test_type", "Type of test to run", NULL, TestTypes, (PetscEnum)(app_ctx->test_type), (PetscEnum *)&app_ctx->test_type, 488fb33541SJames Wright NULL)); 4977841947SLeila Ghaffari 5077841947SLeila Ghaffari app_ctx->test_tol = 1E-11; 512b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-compare_final_state_atol", "Test absolute tolerance", NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL)); 5277841947SLeila Ghaffari 53b7d66439SJames Wright PetscCall(PetscOptionsString("-compare_final_state_filename", "Test filename", NULL, app_ctx->test_file_path, app_ctx->test_file_path, 54b7d66439SJames Wright sizeof(app_ctx->test_file_path), NULL)); 5577841947SLeila Ghaffari 562b730f8bSJeremy L Thompson PetscCall(PetscOptionsFList("-problem", "Problem to solve", NULL, app_ctx->problems, app_ctx->problem_name, app_ctx->problem_name, 572b730f8bSJeremy L Thompson sizeof(app_ctx->problem_name), &problem_flag)); 5877841947SLeila Ghaffari 5977841947SLeila Ghaffari app_ctx->viz_refine = 0; 602b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-viz_refine", "Regular refinement levels for visualization", NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL)); 6177841947SLeila Ghaffari 6237cbb16aSJed Brown app_ctx->checkpoint_interval = 10; 6337cbb16aSJed Brown app_ctx->checkpoint_vtk = PETSC_FALSE; 6437cbb16aSJed Brown PetscCall(PetscOptionsDeprecated("-output_freq", "-checkpoint_interval", "libCEED 0.11.1", "Use -checkpoint_vtk true to include VTK output")); 6537cbb16aSJed Brown PetscCall(PetscOptionsInt("-output_freq", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval, 6637cbb16aSJed Brown &app_ctx->checkpoint_interval, &option_set)); 6737cbb16aSJed Brown if (option_set) app_ctx->checkpoint_vtk = PETSC_TRUE; 6837cbb16aSJed Brown PetscCall(PetscOptionsInt("-checkpoint_interval", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval, 6937cbb16aSJed Brown &app_ctx->checkpoint_interval, NULL)); 7037cbb16aSJed Brown PetscCall(PetscOptionsBool("-checkpoint_vtk", "Include VTK (*.vtu) output at each binary checkpoint", NULL, app_ctx->checkpoint_vtk, 7137cbb16aSJed Brown &app_ctx->checkpoint_vtk, NULL)); 7277841947SLeila Ghaffari 732b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-output_add_stepnum2bin", "Add step number to the binary outputs", NULL, app_ctx->add_stepnum2bin, 742b730f8bSJeremy L Thompson &app_ctx->add_stepnum2bin, NULL)); 7569293791SJames Wright 7669293791SJames Wright PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2)); 772b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL)); 7869293791SJames Wright 7977841947SLeila Ghaffari app_ctx->cont_steps = 0; 802b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-continue", "Continue from previous solution", NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL)); 8177841947SLeila Ghaffari 8269293791SJames Wright PetscCall(PetscStrcpy(app_ctx->cont_file, "[output_dir]/ns-solution.bin")); 832b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-continue_filename", "Filename to get initial condition from", NULL, app_ctx->cont_file, app_ctx->cont_file, 8469293791SJames Wright sizeof(app_ctx->cont_file), &option_set)); 852b730f8bSJeremy L Thompson if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_file, sizeof app_ctx->cont_file, "%s/ns-solution.bin", app_ctx->output_dir)); 864de8550aSJed Brown if (option_set && app_ctx->cont_steps == 0) app_ctx->cont_steps = -1; // Read time from file 8769293791SJames Wright 8869293791SJames Wright PetscCall(PetscStrcpy(app_ctx->cont_time_file, "[output_dir]/ns-time.bin")); 892b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-continue_time_filename", "Filename to get initial condition time from", NULL, app_ctx->cont_time_file, 902b730f8bSJeremy L Thompson app_ctx->cont_time_file, sizeof(app_ctx->cont_time_file), &option_set)); 912b730f8bSJeremy L Thompson if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_time_file, sizeof app_ctx->cont_time_file, "%s/ns-time.bin", app_ctx->output_dir)); 9269293791SJames Wright 9377841947SLeila Ghaffari app_ctx->degree = 1; 942b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of finite elements", NULL, app_ctx->degree, &app_ctx->degree, NULL)); 9577841947SLeila Ghaffari 96fc14f3f6SLeila Ghaffari app_ctx->q_extra = 0; 972b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL)); 9877841947SLeila Ghaffari 99544be873SJed Brown { 100544be873SJed Brown PetscBool option_set; 101544be873SJed Brown char amat_type[256] = ""; 1022b730f8bSJeremy L Thompson PetscCall(PetscOptionsFList("-amat_type", "Set the type of Amat distinct from Pmat (-dm_mat_type)", NULL, MatList, amat_type, amat_type, 1032b730f8bSJeremy L Thompson sizeof(amat_type), &option_set)); 1042b730f8bSJeremy L Thompson if (option_set) PetscCall(PetscStrallocpy(amat_type, (char **)&app_ctx->amat_type)); 105544be873SJed Brown } 1062b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-pmat_pbdiagonal", "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal, 107544be873SJed Brown &app_ctx->pmat_pbdiagonal, NULL)); 108544be873SJed Brown 10977841947SLeila Ghaffari // Provide default ceed resource if not specified 11077841947SLeila Ghaffari if (!ceed_flag) { 11177841947SLeila Ghaffari const char *ceed_resource = "/cpu/self"; 11277841947SLeila Ghaffari strncpy(app_ctx->ceed_resource, ceed_resource, 10); 11377841947SLeila Ghaffari } 1140992c3d1SJed Brown // If we request a GPU, make sure PETSc has initialized its device (which is 1150992c3d1SJed Brown // MPI-aware in case multiple devices are available) before CeedInit so that 1160992c3d1SJed Brown // PETSc and libCEED agree about which device to use. 1170992c3d1SJed Brown if (strncmp(app_ctx->ceed_resource, "/gpu", 4) == 0) PetscCall(PetscDeviceInitialize(PETSC_DEVICE_DEFAULT())); 11877841947SLeila Ghaffari 11977841947SLeila Ghaffari // Provide default problem if not specified 12077841947SLeila Ghaffari if (!problem_flag) { 12177841947SLeila Ghaffari const char *problem_name = "density_current"; 12277841947SLeila Ghaffari strncpy(app_ctx->problem_name, problem_name, 16); 12377841947SLeila Ghaffari } 12477841947SLeila Ghaffari 1252fe7aee7SLeila Ghaffari // Wall Boundary Conditions 1262fe7aee7SLeila Ghaffari bc->num_wall = 16; 1272fe7aee7SLeila Ghaffari PetscBool flg; 1282b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_wall", "Face IDs to apply wall BC", NULL, bc->walls, &bc->num_wall, NULL)); 1292fe7aee7SLeila Ghaffari bc->num_comps = 5; 1302b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-wall_comps", "An array of constrained component numbers", NULL, bc->wall_comps, &bc->num_comps, &flg)); 1312fe7aee7SLeila Ghaffari // Slip Boundary Conditions 1322fe7aee7SLeila Ghaffari for (PetscInt j = 0; j < 3; j++) { 1332fe7aee7SLeila Ghaffari bc->num_slip[j] = 16; 1342fe7aee7SLeila Ghaffari PetscBool flg; 1352fe7aee7SLeila Ghaffari const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1362b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray(flags[j], "Face IDs to apply slip BC", NULL, bc->slips[j], &bc->num_slip[j], &flg)); 1372fe7aee7SLeila Ghaffari if (flg) bc->user_bc = PETSC_TRUE; 1382fe7aee7SLeila Ghaffari } 1392fe7aee7SLeila Ghaffari 1402fe7aee7SLeila Ghaffari // Error if wall and slip BCs are set on the same face 14149aac155SJeremy L Thompson if (bc->user_bc) { 14249aac155SJeremy L Thompson for (PetscInt c = 0; c < 3; c++) { 14349aac155SJeremy L Thompson for (PetscInt s = 0; s < bc->num_slip[c]; s++) { 14449aac155SJeremy L Thompson for (PetscInt w = 0; w < bc->num_wall; w++) { 1450e654f56SJames Wright PetscCheck(bc->slips[c][s] != bc->walls[w], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 1460e654f56SJames Wright "Boundary condition already set on face %" PetscInt_FMT "!\n", bc->walls[w]); 14749aac155SJeremy L Thompson } 14849aac155SJeremy L Thompson } 14949aac155SJeremy L Thompson } 15049aac155SJeremy L Thompson } 151f5452247SJames Wright app_ctx->wall_forces.num_wall = bc->num_wall; 152f5452247SJames Wright PetscMalloc1(bc->num_wall, &app_ctx->wall_forces.walls); 153f5452247SJames Wright PetscCall(PetscArraycpy(app_ctx->wall_forces.walls, bc->walls, bc->num_wall)); 1542fe7aee7SLeila Ghaffari 1552fe7aee7SLeila Ghaffari // Inflow BCs 1562fe7aee7SLeila Ghaffari bc->num_inflow = 16; 1572b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_inflow", "Face IDs to apply inflow BC", NULL, bc->inflows, &bc->num_inflow, NULL)); 1582fe7aee7SLeila Ghaffari // Outflow BCs 1592fe7aee7SLeila Ghaffari bc->num_outflow = 16; 1602b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_outflow", "Face IDs to apply outflow BC", NULL, bc->outflows, &bc->num_outflow, NULL)); 161f4ca79c2SJames Wright // Freestream BCs 162f4ca79c2SJames Wright bc->num_freestream = 16; 1632b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_freestream", "Face IDs to apply freestream BC", NULL, bc->freestreams, &bc->num_freestream, NULL)); 1642fe7aee7SLeila Ghaffari 165a175e481SJames Wright // Statistics Options 166b7d66439SJames Wright app_ctx->turb_spanstats_collect_interval = 1; 167b7d66439SJames Wright PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_collect_interval", "Number of timesteps between statistics collection", NULL, 168b7d66439SJames Wright app_ctx->turb_spanstats_collect_interval, &app_ctx->turb_spanstats_collect_interval, NULL)); 169a175e481SJames Wright 170b7d66439SJames Wright app_ctx->turb_spanstats_viewer_interval = -1; 171b7d66439SJames Wright PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_viewer_interval", "Number of timesteps between statistics viewer writing", NULL, 172b7d66439SJames Wright app_ctx->turb_spanstats_viewer_interval, &app_ctx->turb_spanstats_viewer_interval, NULL)); 173a175e481SJames Wright 174b7d66439SJames Wright PetscCall(PetscOptionsViewer("-ts_monitor_turbulence_spanstats_viewer", "Viewer for the statistics", NULL, &app_ctx->turb_spanstats_viewer, 175b7d66439SJames Wright &app_ctx->turb_spanstats_viewer_format, &app_ctx->turb_spanstats_enable)); 176a175e481SJames Wright 177ca69d878SAdeleke O. Bankole PetscCall(PetscOptionsViewer("-ts_monitor_wall_force", "Viewer for force on each (no-slip) wall", NULL, &app_ctx->wall_forces.viewer, 178ca69d878SAdeleke O. Bankole &app_ctx->wall_forces.viewer_format, NULL)); 179ca69d878SAdeleke O. Bankole 18019ffbc25SJames Wright // SGS Model Options 18119ffbc25SJames Wright app_ctx->sgs_model_type = SGS_MODEL_NONE; 18219ffbc25SJames Wright PetscCall(PetscOptionsEnum("-sgs_model_type", "Subgrid Stress Model type", NULL, SGSModelTypes, (PetscEnum)app_ctx->sgs_model_type, 18319ffbc25SJames Wright (PetscEnum *)&app_ctx->sgs_model_type, NULL)); 18419ffbc25SJames Wright 1854e9802d1SJames Wright PetscCall(PetscOptionsBool("-diff_filter_monitor", "Enable differential filtering TSMonitor", NULL, app_ctx->diff_filter_monitor, 1864e9802d1SJames Wright &app_ctx->diff_filter_monitor, NULL)); 1874e9802d1SJames Wright 1882526956eSJames Wright // Mesh Transformation Options 1892526956eSJames Wright app_ctx->mesh_transform_type = MESH_TRANSFORM_NONE; 1902526956eSJames Wright PetscCall(PetscOptionsEnum("-mesh_transform", "Mesh transform to perform", NULL, MeshTransformTypes, (PetscEnum)app_ctx->mesh_transform_type, 1912526956eSJames Wright (PetscEnum *)&app_ctx->mesh_transform_type, NULL)); 1922526956eSJames Wright 19367490bc6SJeremy L Thompson PetscOptionsEnd(); 194ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 19577841947SLeila Ghaffari } 196