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 1177841947SLeila Ghaffari #include "../navierstokes.h" 1277841947SLeila Ghaffari 1377841947SLeila Ghaffari // Register problems to be available on the command line 1477841947SLeila Ghaffari PetscErrorCode RegisterProblems_NS(AppCtx app_ctx) { 1577841947SLeila Ghaffari 1677841947SLeila Ghaffari 1777841947SLeila Ghaffari app_ctx->problems = NULL; 1877841947SLeila Ghaffari PetscErrorCode ierr; 1977841947SLeila Ghaffari PetscFunctionBeginUser; 2077841947SLeila Ghaffari 2177841947SLeila Ghaffari ierr = PetscFunctionListAdd(&app_ctx->problems, "density_current", 2277841947SLeila Ghaffari NS_DENSITY_CURRENT); CHKERRQ(ierr); 2377841947SLeila Ghaffari 2477841947SLeila Ghaffari ierr = PetscFunctionListAdd(&app_ctx->problems, "euler_vortex", 2577841947SLeila Ghaffari NS_EULER_VORTEX); CHKERRQ(ierr); 2677841947SLeila Ghaffari 27019b7682STimothy Aiken ierr = PetscFunctionListAdd(&app_ctx->problems, "shocktube", 28019b7682STimothy Aiken NS_SHOCKTUBE); CHKERRQ(ierr); 29019b7682STimothy Aiken 3077841947SLeila Ghaffari ierr = PetscFunctionListAdd(&app_ctx->problems, "advection", 3177841947SLeila Ghaffari NS_ADVECTION); CHKERRQ(ierr); 3277841947SLeila Ghaffari 3377841947SLeila Ghaffari ierr = PetscFunctionListAdd(&app_ctx->problems, "advection2d", 3477841947SLeila Ghaffari NS_ADVECTION2D); CHKERRQ(ierr); 3577841947SLeila Ghaffari 3688626eedSJames Wright ierr = PetscFunctionListAdd(&app_ctx->problems, "blasius", 3788626eedSJames Wright NS_BLASIUS); CHKERRQ(ierr); 3888626eedSJames Wright 3988626eedSJames Wright ierr = PetscFunctionListAdd(&app_ctx->problems, "channel", 4088626eedSJames Wright NS_CHANNEL); CHKERRQ(ierr); 4188626eedSJames Wright 4277841947SLeila Ghaffari PetscFunctionReturn(0); 4377841947SLeila Ghaffari } 4477841947SLeila Ghaffari 4577841947SLeila Ghaffari // Process general command line options 462fe7aee7SLeila Ghaffari PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, 472fe7aee7SLeila Ghaffari SimpleBC bc) { 4877841947SLeila Ghaffari 4977841947SLeila Ghaffari PetscBool ceed_flag = PETSC_FALSE; 5077841947SLeila Ghaffari PetscBool problem_flag = PETSC_FALSE; 5177841947SLeila Ghaffari PetscErrorCode ierr; 5277841947SLeila Ghaffari PetscFunctionBeginUser; 5377841947SLeila Ghaffari 5467490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", 5567490bc6SJeremy L Thompson NULL); 5677841947SLeila Ghaffari 5777841947SLeila Ghaffari ierr = PetscOptionsString("-ceed", "CEED resource specifier", 5877841947SLeila Ghaffari NULL, app_ctx->ceed_resource, app_ctx->ceed_resource, 5977841947SLeila Ghaffari sizeof(app_ctx->ceed_resource), &ceed_flag); CHKERRQ(ierr); 6077841947SLeila Ghaffari 6177841947SLeila Ghaffari app_ctx->test_mode = PETSC_FALSE; 6277841947SLeila Ghaffari ierr = PetscOptionsBool("-test", "Run in test mode", 6377841947SLeila Ghaffari NULL, app_ctx->test_mode, &app_ctx->test_mode, NULL); CHKERRQ(ierr); 6477841947SLeila Ghaffari 6577841947SLeila Ghaffari app_ctx->test_tol = 1E-11; 6677841947SLeila Ghaffari ierr = PetscOptionsScalar("-compare_final_state_atol", 6777841947SLeila Ghaffari "Test absolute tolerance", 6877841947SLeila Ghaffari NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL); CHKERRQ(ierr); 6977841947SLeila Ghaffari 7077841947SLeila Ghaffari ierr = PetscOptionsString("-compare_final_state_filename", "Test filename", 7177841947SLeila Ghaffari NULL, app_ctx->file_path, app_ctx->file_path, 7277841947SLeila Ghaffari sizeof(app_ctx->file_path), NULL); CHKERRQ(ierr); 7377841947SLeila Ghaffari 7477841947SLeila Ghaffari ierr = PetscOptionsFList("-problem", "Problem to solve", NULL, 7577841947SLeila Ghaffari app_ctx->problems, 7677841947SLeila Ghaffari app_ctx->problem_name, app_ctx->problem_name, sizeof(app_ctx->problem_name), 7777841947SLeila Ghaffari &problem_flag); CHKERRQ(ierr); 7877841947SLeila Ghaffari 7977841947SLeila Ghaffari app_ctx->viz_refine = 0; 8077841947SLeila Ghaffari ierr = PetscOptionsInt("-viz_refine", 8177841947SLeila Ghaffari "Regular refinement levels for visualization", 8277841947SLeila Ghaffari NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL); CHKERRQ(ierr); 8377841947SLeila Ghaffari 8477841947SLeila Ghaffari app_ctx->output_freq = 10; 8577841947SLeila Ghaffari ierr = PetscOptionsInt("-output_freq", 8677841947SLeila Ghaffari "Frequency of output, in number of steps", 8777841947SLeila Ghaffari NULL, app_ctx->output_freq, &app_ctx->output_freq, NULL); CHKERRQ(ierr); 8877841947SLeila Ghaffari 8977841947SLeila Ghaffari app_ctx->cont_steps = 0; 9077841947SLeila Ghaffari ierr = PetscOptionsInt("-continue", "Continue from previous solution", 9177841947SLeila Ghaffari NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL); CHKERRQ(ierr); 9277841947SLeila Ghaffari 9377841947SLeila Ghaffari app_ctx->degree = 1; 9477841947SLeila Ghaffari ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements", 9577841947SLeila Ghaffari NULL, app_ctx->degree, &app_ctx->degree, NULL); CHKERRQ(ierr); 9677841947SLeila Ghaffari 9777841947SLeila Ghaffari app_ctx->q_extra = 2; 9877841947SLeila Ghaffari ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points", 9977841947SLeila Ghaffari NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL); CHKERRQ(ierr); 10077841947SLeila Ghaffari 101544be873SJed Brown { 102544be873SJed Brown PetscBool option_set; 103544be873SJed Brown char amat_type[256] = ""; 104544be873SJed Brown PetscCall(PetscOptionsFList("-amat_type", 105544be873SJed Brown "Set the type of Amat distinct from Pmat (-dm_mat_type)", 106544be873SJed Brown NULL, MatList, amat_type, amat_type, sizeof(amat_type), &option_set)); 107544be873SJed Brown if (option_set) PetscCall(PetscStrallocpy(amat_type, 108544be873SJed Brown (char **)&app_ctx->amat_type)); 109544be873SJed Brown } 110544be873SJed Brown PetscCall(PetscOptionsBool("-pmat_pbdiagonal", 111544be873SJed Brown "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal, 112544be873SJed Brown &app_ctx->pmat_pbdiagonal, NULL)); 113544be873SJed Brown 11477841947SLeila Ghaffari ierr = PetscStrncpy(app_ctx->output_dir, ".", 2); CHKERRQ(ierr); 11577841947SLeila Ghaffari ierr = PetscOptionsString("-output_dir", "Output directory", 11677841947SLeila Ghaffari NULL, app_ctx->output_dir, app_ctx->output_dir, 11777841947SLeila Ghaffari sizeof(app_ctx->output_dir), NULL); CHKERRQ(ierr); 11877841947SLeila Ghaffari 11977841947SLeila Ghaffari // Provide default ceed resource if not specified 12077841947SLeila Ghaffari if (!ceed_flag) { 12177841947SLeila Ghaffari const char *ceed_resource = "/cpu/self"; 12277841947SLeila Ghaffari strncpy(app_ctx->ceed_resource, ceed_resource, 10); 12377841947SLeila Ghaffari } 12477841947SLeila Ghaffari 12577841947SLeila Ghaffari // Provide default problem if not specified 12677841947SLeila Ghaffari if (!problem_flag) { 12777841947SLeila Ghaffari const char *problem_name = "density_current"; 12877841947SLeila Ghaffari strncpy(app_ctx->problem_name, problem_name, 16); 12977841947SLeila Ghaffari } 13077841947SLeila Ghaffari 1312fe7aee7SLeila Ghaffari // Wall Boundary Conditions 1322fe7aee7SLeila Ghaffari bc->num_wall = 16; 1332fe7aee7SLeila Ghaffari PetscBool flg; 1342fe7aee7SLeila Ghaffari ierr = PetscOptionsIntArray("-bc_wall", 1352fe7aee7SLeila Ghaffari "Face IDs to apply wall BC", 1362fe7aee7SLeila Ghaffari NULL, bc->walls, &bc->num_wall, NULL); CHKERRQ(ierr); 1372fe7aee7SLeila Ghaffari bc->num_comps = 5; 1382fe7aee7SLeila Ghaffari ierr = PetscOptionsIntArray("-wall_comps", 1392fe7aee7SLeila Ghaffari "An array of constrained component numbers", 1402fe7aee7SLeila Ghaffari NULL, bc->wall_comps, &bc->num_comps, &flg); CHKERRQ(ierr); 1412fe7aee7SLeila Ghaffari // Slip Boundary Conditions 1422fe7aee7SLeila Ghaffari for (PetscInt j=0; j<3; j++) { 1432fe7aee7SLeila Ghaffari bc->num_slip[j] = 16; 1442fe7aee7SLeila Ghaffari PetscBool flg; 1452fe7aee7SLeila Ghaffari const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"}; 1462fe7aee7SLeila Ghaffari ierr = PetscOptionsIntArray(flags[j], 1472fe7aee7SLeila Ghaffari "Face IDs to apply slip BC", 1482fe7aee7SLeila Ghaffari NULL, bc->slips[j], &bc->num_slip[j], &flg); CHKERRQ(ierr); 1492fe7aee7SLeila Ghaffari if (flg) bc->user_bc = PETSC_TRUE; 1502fe7aee7SLeila Ghaffari } 1512fe7aee7SLeila Ghaffari 1522fe7aee7SLeila Ghaffari // Error if wall and slip BCs are set on the same face 1532fe7aee7SLeila Ghaffari if (bc->user_bc) 1542fe7aee7SLeila Ghaffari for (PetscInt c = 0; c < 3; c++) 1552fe7aee7SLeila Ghaffari for (PetscInt s = 0; s < bc->num_slip[c]; s++) 1562fe7aee7SLeila Ghaffari for (PetscInt w = 0; w < bc->num_wall; w++) 1572fe7aee7SLeila Ghaffari if (bc->slips[c][s] == bc->walls[w]) 1587578c821SJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, 15908140895SJed Brown "Boundary condition already set on face %" PetscInt_FMT "!\n", 1602fe7aee7SLeila Ghaffari bc->walls[w]); 1612fe7aee7SLeila Ghaffari 1622fe7aee7SLeila Ghaffari // Inflow BCs 1632fe7aee7SLeila Ghaffari bc->num_inflow = 16; 1642fe7aee7SLeila Ghaffari ierr = PetscOptionsIntArray("-bc_inflow", 1652fe7aee7SLeila Ghaffari "Face IDs to apply inflow BC", 1662fe7aee7SLeila Ghaffari NULL, bc->inflows, &bc->num_inflow, NULL); CHKERRQ(ierr); 1672fe7aee7SLeila Ghaffari // Outflow BCs 1682fe7aee7SLeila Ghaffari bc->num_outflow = 16; 1692fe7aee7SLeila Ghaffari ierr = PetscOptionsIntArray("-bc_outflow", 1702fe7aee7SLeila Ghaffari "Face IDs to apply outflow BC", 1712fe7aee7SLeila Ghaffari NULL, bc->outflows, &bc->num_outflow, NULL); CHKERRQ(ierr); 172*f4ca79c2SJames Wright // Freestream BCs 173*f4ca79c2SJames Wright bc->num_freestream = 16; 174*f4ca79c2SJames Wright ierr = PetscOptionsIntArray("-bc_freestream", 175*f4ca79c2SJames Wright "Face IDs to apply freestream BC", 176*f4ca79c2SJames Wright NULL, bc->freestreams, &bc->num_freestream, NULL); CHKERRQ(ierr); 1772fe7aee7SLeila Ghaffari 17867490bc6SJeremy L Thompson PetscOptionsEnd(); 17977841947SLeila Ghaffari 18077841947SLeila Ghaffari PetscFunctionReturn(0); 18177841947SLeila Ghaffari } 182