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> 12*49aac155SJeremy 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 PetscFunctionBeginUser; 2077841947SLeila Ghaffari 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)); 3088626eedSJames Wright 3177841947SLeila Ghaffari PetscFunctionReturn(0); 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; 4177841947SLeila Ghaffari 422b730f8bSJeremy L Thompson PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", NULL); 4377841947SLeila Ghaffari 442b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource, 452b730f8bSJeremy L Thompson sizeof(app_ctx->ceed_resource), &ceed_flag)); 4677841947SLeila Ghaffari 478fb33541SJames Wright app_ctx->test_type = TESTTYPE_NONE; 488fb33541SJames Wright PetscCall(PetscOptionsEnum("-test_type", "Type of test to run", NULL, TestTypes, (PetscEnum)(app_ctx->test_type), (PetscEnum *)&app_ctx->test_type, 498fb33541SJames Wright NULL)); 5077841947SLeila Ghaffari 5177841947SLeila Ghaffari app_ctx->test_tol = 1E-11; 522b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-compare_final_state_atol", "Test absolute tolerance", NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL)); 5377841947SLeila Ghaffari 54b7d66439SJames Wright PetscCall(PetscOptionsString("-compare_final_state_filename", "Test filename", NULL, app_ctx->test_file_path, app_ctx->test_file_path, 55b7d66439SJames Wright sizeof(app_ctx->test_file_path), NULL)); 5677841947SLeila Ghaffari 572b730f8bSJeremy L Thompson PetscCall(PetscOptionsFList("-problem", "Problem to solve", NULL, app_ctx->problems, app_ctx->problem_name, app_ctx->problem_name, 582b730f8bSJeremy L Thompson sizeof(app_ctx->problem_name), &problem_flag)); 5977841947SLeila Ghaffari 6077841947SLeila Ghaffari app_ctx->viz_refine = 0; 612b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-viz_refine", "Regular refinement levels for visualization", NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL)); 6277841947SLeila Ghaffari 6337cbb16aSJed Brown app_ctx->checkpoint_interval = 10; 6437cbb16aSJed Brown app_ctx->checkpoint_vtk = PETSC_FALSE; 6537cbb16aSJed Brown PetscCall(PetscOptionsDeprecated("-output_freq", "-checkpoint_interval", "libCEED 0.11.1", "Use -checkpoint_vtk true to include VTK output")); 6637cbb16aSJed Brown PetscCall(PetscOptionsInt("-output_freq", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval, 6737cbb16aSJed Brown &app_ctx->checkpoint_interval, &option_set)); 6837cbb16aSJed Brown if (option_set) app_ctx->checkpoint_vtk = PETSC_TRUE; 6937cbb16aSJed Brown PetscCall(PetscOptionsInt("-checkpoint_interval", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval, 7037cbb16aSJed Brown &app_ctx->checkpoint_interval, NULL)); 7137cbb16aSJed Brown PetscCall(PetscOptionsBool("-checkpoint_vtk", "Include VTK (*.vtu) output at each binary checkpoint", NULL, app_ctx->checkpoint_vtk, 7237cbb16aSJed Brown &app_ctx->checkpoint_vtk, NULL)); 7377841947SLeila Ghaffari 742b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-output_add_stepnum2bin", "Add step number to the binary outputs", NULL, app_ctx->add_stepnum2bin, 752b730f8bSJeremy L Thompson &app_ctx->add_stepnum2bin, NULL)); 7669293791SJames Wright 7769293791SJames Wright PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2)); 782b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL)); 7969293791SJames Wright 8077841947SLeila Ghaffari app_ctx->cont_steps = 0; 812b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-continue", "Continue from previous solution", NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL)); 8277841947SLeila Ghaffari 8369293791SJames Wright PetscCall(PetscStrcpy(app_ctx->cont_file, "[output_dir]/ns-solution.bin")); 842b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-continue_filename", "Filename to get initial condition from", NULL, app_ctx->cont_file, app_ctx->cont_file, 8569293791SJames Wright sizeof(app_ctx->cont_file), &option_set)); 862b730f8bSJeremy L Thompson if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_file, sizeof app_ctx->cont_file, "%s/ns-solution.bin", app_ctx->output_dir)); 874de8550aSJed Brown if (option_set && app_ctx->cont_steps == 0) app_ctx->cont_steps = -1; // Read time from file 8869293791SJames Wright 8969293791SJames Wright PetscCall(PetscStrcpy(app_ctx->cont_time_file, "[output_dir]/ns-time.bin")); 902b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-continue_time_filename", "Filename to get initial condition time from", NULL, app_ctx->cont_time_file, 912b730f8bSJeremy L Thompson app_ctx->cont_time_file, sizeof(app_ctx->cont_time_file), &option_set)); 922b730f8bSJeremy 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)); 9369293791SJames Wright 9477841947SLeila Ghaffari app_ctx->degree = 1; 952b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of finite elements", NULL, app_ctx->degree, &app_ctx->degree, NULL)); 9677841947SLeila Ghaffari 97fc14f3f6SLeila Ghaffari app_ctx->q_extra = 0; 982b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL)); 9977841947SLeila Ghaffari 100544be873SJed Brown { 101544be873SJed Brown PetscBool option_set; 102544be873SJed Brown char amat_type[256] = ""; 1032b730f8bSJeremy L Thompson PetscCall(PetscOptionsFList("-amat_type", "Set the type of Amat distinct from Pmat (-dm_mat_type)", NULL, MatList, amat_type, amat_type, 1042b730f8bSJeremy L Thompson sizeof(amat_type), &option_set)); 1052b730f8bSJeremy L Thompson if (option_set) PetscCall(PetscStrallocpy(amat_type, (char **)&app_ctx->amat_type)); 106544be873SJed Brown } 1072b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-pmat_pbdiagonal", "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal, 108544be873SJed Brown &app_ctx->pmat_pbdiagonal, NULL)); 109544be873SJed Brown 11077841947SLeila Ghaffari // Provide default ceed resource if not specified 11177841947SLeila Ghaffari if (!ceed_flag) { 11277841947SLeila Ghaffari const char *ceed_resource = "/cpu/self"; 11377841947SLeila Ghaffari strncpy(app_ctx->ceed_resource, ceed_resource, 10); 11477841947SLeila Ghaffari } 1150992c3d1SJed Brown // If we request a GPU, make sure PETSc has initialized its device (which is 1160992c3d1SJed Brown // MPI-aware in case multiple devices are available) before CeedInit so that 1170992c3d1SJed Brown // PETSc and libCEED agree about which device to use. 1180992c3d1SJed Brown if (strncmp(app_ctx->ceed_resource, "/gpu", 4) == 0) PetscCall(PetscDeviceInitialize(PETSC_DEVICE_DEFAULT())); 11977841947SLeila Ghaffari 12077841947SLeila Ghaffari // Provide default problem if not specified 12177841947SLeila Ghaffari if (!problem_flag) { 12277841947SLeila Ghaffari const char *problem_name = "density_current"; 12377841947SLeila Ghaffari strncpy(app_ctx->problem_name, problem_name, 16); 12477841947SLeila Ghaffari } 12577841947SLeila Ghaffari 1262fe7aee7SLeila Ghaffari // Wall Boundary Conditions 1272fe7aee7SLeila Ghaffari bc->num_wall = 16; 1282fe7aee7SLeila Ghaffari PetscBool flg; 1292b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_wall", "Face IDs to apply wall BC", NULL, bc->walls, &bc->num_wall, NULL)); 1302fe7aee7SLeila Ghaffari bc->num_comps = 5; 1312b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-wall_comps", "An array of constrained component numbers", NULL, bc->wall_comps, &bc->num_comps, &flg)); 1322fe7aee7SLeila Ghaffari // Slip Boundary Conditions 1332fe7aee7SLeila Ghaffari for (PetscInt j = 0; j < 3; j++) { 1342fe7aee7SLeila Ghaffari bc->num_slip[j] = 16; 1352fe7aee7SLeila Ghaffari PetscBool flg; 1362fe7aee7SLeila Ghaffari const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1372b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray(flags[j], "Face IDs to apply slip BC", NULL, bc->slips[j], &bc->num_slip[j], &flg)); 1382fe7aee7SLeila Ghaffari if (flg) bc->user_bc = PETSC_TRUE; 1392fe7aee7SLeila Ghaffari } 1402fe7aee7SLeila Ghaffari 1412fe7aee7SLeila Ghaffari // Error if wall and slip BCs are set on the same face 142*49aac155SJeremy L Thompson if (bc->user_bc) { 143*49aac155SJeremy L Thompson for (PetscInt c = 0; c < 3; c++) { 144*49aac155SJeremy L Thompson for (PetscInt s = 0; s < bc->num_slip[c]; s++) { 145*49aac155SJeremy L Thompson for (PetscInt w = 0; w < bc->num_wall; w++) { 146*49aac155SJeremy L Thompson if (bc->slips[c][s] == bc->walls[w]) { 1472b730f8bSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary condition already set on face %" PetscInt_FMT "!\n", bc->walls[w]); 148*49aac155SJeremy L Thompson } 149*49aac155SJeremy L Thompson } 150*49aac155SJeremy L Thompson } 151*49aac155SJeremy L Thompson } 152*49aac155SJeremy L Thompson } 1532fe7aee7SLeila Ghaffari 1542fe7aee7SLeila Ghaffari // Inflow BCs 1552fe7aee7SLeila Ghaffari bc->num_inflow = 16; 1562b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_inflow", "Face IDs to apply inflow BC", NULL, bc->inflows, &bc->num_inflow, NULL)); 1572fe7aee7SLeila Ghaffari // Outflow BCs 1582fe7aee7SLeila Ghaffari bc->num_outflow = 16; 1592b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_outflow", "Face IDs to apply outflow BC", NULL, bc->outflows, &bc->num_outflow, NULL)); 160f4ca79c2SJames Wright // Freestream BCs 161f4ca79c2SJames Wright bc->num_freestream = 16; 1622b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-bc_freestream", "Face IDs to apply freestream BC", NULL, bc->freestreams, &bc->num_freestream, NULL)); 1632fe7aee7SLeila Ghaffari 164a175e481SJames Wright // Statistics Options 165b7d66439SJames Wright app_ctx->turb_spanstats_collect_interval = 1; 166b7d66439SJames Wright PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_collect_interval", "Number of timesteps between statistics collection", NULL, 167b7d66439SJames Wright app_ctx->turb_spanstats_collect_interval, &app_ctx->turb_spanstats_collect_interval, NULL)); 168a175e481SJames Wright 169b7d66439SJames Wright app_ctx->turb_spanstats_viewer_interval = -1; 170b7d66439SJames Wright PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_viewer_interval", "Number of timesteps between statistics viewer writing", NULL, 171b7d66439SJames Wright app_ctx->turb_spanstats_viewer_interval, &app_ctx->turb_spanstats_viewer_interval, NULL)); 172a175e481SJames Wright 173b7d66439SJames Wright PetscCall(PetscOptionsViewer("-ts_monitor_turbulence_spanstats_viewer", "Viewer for the statistics", NULL, &app_ctx->turb_spanstats_viewer, 174b7d66439SJames Wright &app_ctx->turb_spanstats_viewer_format, &app_ctx->turb_spanstats_enable)); 175a175e481SJames Wright 176ca69d878SAdeleke O. Bankole PetscCall(PetscOptionsViewer("-ts_monitor_wall_force", "Viewer for force on each (no-slip) wall", NULL, &app_ctx->wall_forces.viewer, 177ca69d878SAdeleke O. Bankole &app_ctx->wall_forces.viewer_format, NULL)); 178ca69d878SAdeleke O. Bankole 17967490bc6SJeremy L Thompson PetscOptionsEnd(); 18077841947SLeila Ghaffari 18177841947SLeila Ghaffari PetscFunctionReturn(0); 18277841947SLeila Ghaffari } 183