xref: /libCEED/examples/fluids/src/cloptions.c (revision 67490bc6de48fa3615ca135c74fbc2a7c9584dee)
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 
2777841947SLeila Ghaffari   ierr = PetscFunctionListAdd(&app_ctx->problems, "advection",
2877841947SLeila Ghaffari                               NS_ADVECTION); CHKERRQ(ierr);
2977841947SLeila Ghaffari 
3077841947SLeila Ghaffari   ierr = PetscFunctionListAdd(&app_ctx->problems, "advection2d",
3177841947SLeila Ghaffari                               NS_ADVECTION2D); CHKERRQ(ierr);
3277841947SLeila Ghaffari 
3377841947SLeila Ghaffari   PetscFunctionReturn(0);
3477841947SLeila Ghaffari }
3577841947SLeila Ghaffari 
3677841947SLeila Ghaffari // Process general command line options
372fe7aee7SLeila Ghaffari PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx,
382fe7aee7SLeila Ghaffari     SimpleBC bc) {
3977841947SLeila Ghaffari 
4077841947SLeila Ghaffari   PetscBool ceed_flag = PETSC_FALSE;
4177841947SLeila Ghaffari   PetscBool problem_flag = PETSC_FALSE;
4277841947SLeila Ghaffari   PetscErrorCode ierr;
4377841947SLeila Ghaffari   PetscFunctionBeginUser;
4477841947SLeila Ghaffari 
45*67490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
46*67490bc6SJeremy L Thompson                     NULL);
4777841947SLeila Ghaffari 
4877841947SLeila Ghaffari   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
4977841947SLeila Ghaffari                             NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
5077841947SLeila Ghaffari                             sizeof(app_ctx->ceed_resource), &ceed_flag); CHKERRQ(ierr);
5177841947SLeila Ghaffari 
5277841947SLeila Ghaffari   app_ctx->test_mode = PETSC_FALSE;
5377841947SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
5477841947SLeila Ghaffari                           NULL, app_ctx->test_mode, &app_ctx->test_mode, NULL); CHKERRQ(ierr);
5577841947SLeila Ghaffari 
5677841947SLeila Ghaffari   app_ctx->test_tol = 1E-11;
5777841947SLeila Ghaffari   ierr = PetscOptionsScalar("-compare_final_state_atol",
5877841947SLeila Ghaffari                             "Test absolute tolerance",
5977841947SLeila Ghaffari                             NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL); CHKERRQ(ierr);
6077841947SLeila Ghaffari 
6177841947SLeila Ghaffari   ierr = PetscOptionsString("-compare_final_state_filename", "Test filename",
6277841947SLeila Ghaffari                             NULL, app_ctx->file_path, app_ctx->file_path,
6377841947SLeila Ghaffari                             sizeof(app_ctx->file_path), NULL); CHKERRQ(ierr);
6477841947SLeila Ghaffari 
6577841947SLeila Ghaffari   ierr = PetscOptionsFList("-problem", "Problem to solve", NULL,
6677841947SLeila Ghaffari                            app_ctx->problems,
6777841947SLeila Ghaffari                            app_ctx->problem_name, app_ctx->problem_name, sizeof(app_ctx->problem_name),
6877841947SLeila Ghaffari                            &problem_flag); CHKERRQ(ierr);
6977841947SLeila Ghaffari 
7077841947SLeila Ghaffari   app_ctx->viz_refine = 0;
7177841947SLeila Ghaffari   ierr = PetscOptionsInt("-viz_refine",
7277841947SLeila Ghaffari                          "Regular refinement levels for visualization",
7377841947SLeila Ghaffari                          NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL); CHKERRQ(ierr);
7477841947SLeila Ghaffari 
7577841947SLeila Ghaffari   app_ctx->output_freq = 10;
7677841947SLeila Ghaffari   ierr = PetscOptionsInt("-output_freq",
7777841947SLeila Ghaffari                          "Frequency of output, in number of steps",
7877841947SLeila Ghaffari                          NULL, app_ctx->output_freq, &app_ctx->output_freq, NULL); CHKERRQ(ierr);
7977841947SLeila Ghaffari 
8077841947SLeila Ghaffari   app_ctx->cont_steps = 0;
8177841947SLeila Ghaffari   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
8277841947SLeila Ghaffari                          NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL); CHKERRQ(ierr);
8377841947SLeila Ghaffari 
8477841947SLeila Ghaffari   app_ctx->degree = 1;
8577841947SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
8677841947SLeila Ghaffari                          NULL, app_ctx->degree, &app_ctx->degree, NULL); CHKERRQ(ierr);
8777841947SLeila Ghaffari 
8877841947SLeila Ghaffari   app_ctx->q_extra = 2;
8977841947SLeila Ghaffari   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
9077841947SLeila Ghaffari                          NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL); CHKERRQ(ierr);
9177841947SLeila Ghaffari 
9277841947SLeila Ghaffari   ierr = PetscStrncpy(app_ctx->output_dir, ".", 2); CHKERRQ(ierr);
9377841947SLeila Ghaffari   ierr = PetscOptionsString("-output_dir", "Output directory",
9477841947SLeila Ghaffari                             NULL, app_ctx->output_dir, app_ctx->output_dir,
9577841947SLeila Ghaffari                             sizeof(app_ctx->output_dir), NULL); CHKERRQ(ierr);
9677841947SLeila Ghaffari 
9777841947SLeila Ghaffari   // Provide default ceed resource if not specified
9877841947SLeila Ghaffari   if (!ceed_flag) {
9977841947SLeila Ghaffari     const char *ceed_resource = "/cpu/self";
10077841947SLeila Ghaffari     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
10177841947SLeila Ghaffari   }
10277841947SLeila Ghaffari 
10377841947SLeila Ghaffari   // Provide default problem if not specified
10477841947SLeila Ghaffari   if (!problem_flag) {
10577841947SLeila Ghaffari     const char *problem_name = "density_current";
10677841947SLeila Ghaffari     strncpy(app_ctx->problem_name, problem_name, 16);
10777841947SLeila Ghaffari   }
10877841947SLeila Ghaffari 
1092fe7aee7SLeila Ghaffari   // Wall Boundary Conditions
1102fe7aee7SLeila Ghaffari   bc->num_wall = 16;
1112fe7aee7SLeila Ghaffari   PetscBool flg;
1122fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_wall",
1132fe7aee7SLeila Ghaffari                               "Face IDs to apply wall BC",
1142fe7aee7SLeila Ghaffari                               NULL, bc->walls, &bc->num_wall, NULL); CHKERRQ(ierr);
1152fe7aee7SLeila Ghaffari   bc->num_comps = 5;
1162fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-wall_comps",
1172fe7aee7SLeila Ghaffari                               "An array of constrained component numbers",
1182fe7aee7SLeila Ghaffari                               NULL, bc->wall_comps, &bc->num_comps, &flg); CHKERRQ(ierr);
1192fe7aee7SLeila Ghaffari   // Slip Boundary Conditions
1202fe7aee7SLeila Ghaffari   for (PetscInt j=0; j<3; j++) {
1212fe7aee7SLeila Ghaffari     bc->num_slip[j] = 16;
1222fe7aee7SLeila Ghaffari     PetscBool flg;
1232fe7aee7SLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1242fe7aee7SLeila Ghaffari     ierr = PetscOptionsIntArray(flags[j],
1252fe7aee7SLeila Ghaffari                                 "Face IDs to apply slip BC",
1262fe7aee7SLeila Ghaffari                                 NULL, bc->slips[j], &bc->num_slip[j], &flg); CHKERRQ(ierr);
1272fe7aee7SLeila Ghaffari     if (flg) bc->user_bc = PETSC_TRUE;
1282fe7aee7SLeila Ghaffari   }
1292fe7aee7SLeila Ghaffari 
1302fe7aee7SLeila Ghaffari   // Error if wall and slip BCs are set on the same face
1312fe7aee7SLeila Ghaffari   if (bc->user_bc)
1322fe7aee7SLeila Ghaffari     for (PetscInt c = 0; c < 3; c++)
1332fe7aee7SLeila Ghaffari       for (PetscInt s = 0; s < bc->num_slip[c]; s++)
1342fe7aee7SLeila Ghaffari         for (PetscInt w = 0; w < bc->num_wall; w++)
1352fe7aee7SLeila Ghaffari           if (bc->slips[c][s] == bc->walls[w])
1367578c821SJeremy L Thompson             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
1372fe7aee7SLeila Ghaffari                     "Boundary condition already set on face %D!\n",
1382fe7aee7SLeila Ghaffari                     bc->walls[w]);
1392fe7aee7SLeila Ghaffari 
1402fe7aee7SLeila Ghaffari   // Inflow BCs
1412fe7aee7SLeila Ghaffari   bc->num_inflow = 16;
1422fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_inflow",
1432fe7aee7SLeila Ghaffari                               "Face IDs to apply inflow BC",
1442fe7aee7SLeila Ghaffari                               NULL, bc->inflows, &bc->num_inflow, NULL); CHKERRQ(ierr);
1452fe7aee7SLeila Ghaffari   // Outflow BCs
1462fe7aee7SLeila Ghaffari   bc->num_outflow = 16;
1472fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_outflow",
1482fe7aee7SLeila Ghaffari                               "Face IDs to apply outflow BC",
1492fe7aee7SLeila Ghaffari                               NULL, bc->outflows, &bc->num_outflow, NULL); CHKERRQ(ierr);
1502fe7aee7SLeila Ghaffari 
151*67490bc6SJeremy L Thompson   PetscOptionsEnd();
15277841947SLeila Ghaffari 
15377841947SLeila Ghaffari   PetscFunctionReturn(0);
15477841947SLeila Ghaffari }
155