xref: /libCEED/examples/fluids/src/cloptions.c (revision cea0a271c0c218b226f18eeac3ff0d23cf73bc2b)
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   PetscFunctionBeginUser;
1977841947SLeila Ghaffari 
20*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "density_current",
21*cea0a271SJames Wright                                  NS_DENSITY_CURRENT));
2277841947SLeila Ghaffari 
23*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "euler_vortex",
24*cea0a271SJames Wright                                  NS_EULER_VORTEX));
2577841947SLeila Ghaffari 
26*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "shocktube",
27*cea0a271SJames Wright                                  NS_SHOCKTUBE));
28019b7682STimothy Aiken 
29*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection",
30*cea0a271SJames Wright                                  NS_ADVECTION));
3177841947SLeila Ghaffari 
32*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection2d",
33*cea0a271SJames Wright                                  NS_ADVECTION2D));
3477841947SLeila Ghaffari 
35*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "blasius",
36*cea0a271SJames Wright                                  NS_BLASIUS));
3788626eedSJames Wright 
38*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "channel",
39*cea0a271SJames Wright                                  NS_CHANNEL));
40*cea0a271SJames Wright 
41*cea0a271SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "newtonian_wave",
42*cea0a271SJames Wright                                  NS_NEWTONIAN_WAVE));
4388626eedSJames Wright 
4477841947SLeila Ghaffari   PetscFunctionReturn(0);
4577841947SLeila Ghaffari }
4677841947SLeila Ghaffari 
4777841947SLeila Ghaffari // Process general command line options
482fe7aee7SLeila Ghaffari PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx,
492fe7aee7SLeila Ghaffari     SimpleBC bc) {
5077841947SLeila Ghaffari 
5177841947SLeila Ghaffari   PetscBool ceed_flag = PETSC_FALSE;
5277841947SLeila Ghaffari   PetscBool problem_flag = PETSC_FALSE;
5377841947SLeila Ghaffari   PetscErrorCode ierr;
5477841947SLeila Ghaffari   PetscFunctionBeginUser;
5577841947SLeila Ghaffari 
5667490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED",
5767490bc6SJeremy L Thompson                     NULL);
5877841947SLeila Ghaffari 
5977841947SLeila Ghaffari   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
6077841947SLeila Ghaffari                             NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
6177841947SLeila Ghaffari                             sizeof(app_ctx->ceed_resource), &ceed_flag); CHKERRQ(ierr);
6277841947SLeila Ghaffari 
6377841947SLeila Ghaffari   app_ctx->test_mode = PETSC_FALSE;
6477841947SLeila Ghaffari   ierr = PetscOptionsBool("-test", "Run in test mode",
6577841947SLeila Ghaffari                           NULL, app_ctx->test_mode, &app_ctx->test_mode, NULL); CHKERRQ(ierr);
6677841947SLeila Ghaffari 
6777841947SLeila Ghaffari   app_ctx->test_tol = 1E-11;
6877841947SLeila Ghaffari   ierr = PetscOptionsScalar("-compare_final_state_atol",
6977841947SLeila Ghaffari                             "Test absolute tolerance",
7077841947SLeila Ghaffari                             NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL); CHKERRQ(ierr);
7177841947SLeila Ghaffari 
7277841947SLeila Ghaffari   ierr = PetscOptionsString("-compare_final_state_filename", "Test filename",
7377841947SLeila Ghaffari                             NULL, app_ctx->file_path, app_ctx->file_path,
7477841947SLeila Ghaffari                             sizeof(app_ctx->file_path), NULL); CHKERRQ(ierr);
7577841947SLeila Ghaffari 
7677841947SLeila Ghaffari   ierr = PetscOptionsFList("-problem", "Problem to solve", NULL,
7777841947SLeila Ghaffari                            app_ctx->problems,
7877841947SLeila Ghaffari                            app_ctx->problem_name, app_ctx->problem_name, sizeof(app_ctx->problem_name),
7977841947SLeila Ghaffari                            &problem_flag); CHKERRQ(ierr);
8077841947SLeila Ghaffari 
8177841947SLeila Ghaffari   app_ctx->viz_refine = 0;
8277841947SLeila Ghaffari   ierr = PetscOptionsInt("-viz_refine",
8377841947SLeila Ghaffari                          "Regular refinement levels for visualization",
8477841947SLeila Ghaffari                          NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL); CHKERRQ(ierr);
8577841947SLeila Ghaffari 
8677841947SLeila Ghaffari   app_ctx->output_freq = 10;
8777841947SLeila Ghaffari   ierr = PetscOptionsInt("-output_freq",
8877841947SLeila Ghaffari                          "Frequency of output, in number of steps",
8977841947SLeila Ghaffari                          NULL, app_ctx->output_freq, &app_ctx->output_freq, NULL); CHKERRQ(ierr);
9077841947SLeila Ghaffari 
9177841947SLeila Ghaffari   app_ctx->cont_steps = 0;
9277841947SLeila Ghaffari   ierr = PetscOptionsInt("-continue", "Continue from previous solution",
9377841947SLeila Ghaffari                          NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL); CHKERRQ(ierr);
9477841947SLeila Ghaffari 
9577841947SLeila Ghaffari   app_ctx->degree = 1;
9677841947SLeila Ghaffari   ierr = PetscOptionsInt("-degree", "Polynomial degree of finite elements",
9777841947SLeila Ghaffari                          NULL, app_ctx->degree, &app_ctx->degree, NULL); CHKERRQ(ierr);
9877841947SLeila Ghaffari 
9977841947SLeila Ghaffari   app_ctx->q_extra = 2;
10077841947SLeila Ghaffari   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
10177841947SLeila Ghaffari                          NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL); CHKERRQ(ierr);
10277841947SLeila Ghaffari 
103544be873SJed Brown   {
104544be873SJed Brown     PetscBool option_set;
105544be873SJed Brown     char amat_type[256] = "";
106544be873SJed Brown     PetscCall(PetscOptionsFList("-amat_type",
107544be873SJed Brown                                 "Set the type of Amat distinct from Pmat (-dm_mat_type)",
108544be873SJed Brown                                 NULL, MatList, amat_type, amat_type, sizeof(amat_type), &option_set));
109544be873SJed Brown     if (option_set) PetscCall(PetscStrallocpy(amat_type,
110544be873SJed Brown                                 (char **)&app_ctx->amat_type));
111544be873SJed Brown   }
112544be873SJed Brown   PetscCall(PetscOptionsBool("-pmat_pbdiagonal",
113544be873SJed Brown                              "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal,
114544be873SJed Brown                              &app_ctx->pmat_pbdiagonal, NULL));
115544be873SJed Brown 
11677841947SLeila Ghaffari   ierr = PetscStrncpy(app_ctx->output_dir, ".", 2); CHKERRQ(ierr);
11777841947SLeila Ghaffari   ierr = PetscOptionsString("-output_dir", "Output directory",
11877841947SLeila Ghaffari                             NULL, app_ctx->output_dir, app_ctx->output_dir,
11977841947SLeila Ghaffari                             sizeof(app_ctx->output_dir), NULL); CHKERRQ(ierr);
12077841947SLeila Ghaffari 
12177841947SLeila Ghaffari   // Provide default ceed resource if not specified
12277841947SLeila Ghaffari   if (!ceed_flag) {
12377841947SLeila Ghaffari     const char *ceed_resource = "/cpu/self";
12477841947SLeila Ghaffari     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
12577841947SLeila Ghaffari   }
12677841947SLeila Ghaffari 
12777841947SLeila Ghaffari   // Provide default problem if not specified
12877841947SLeila Ghaffari   if (!problem_flag) {
12977841947SLeila Ghaffari     const char *problem_name = "density_current";
13077841947SLeila Ghaffari     strncpy(app_ctx->problem_name, problem_name, 16);
13177841947SLeila Ghaffari   }
13277841947SLeila Ghaffari 
1332fe7aee7SLeila Ghaffari   // Wall Boundary Conditions
1342fe7aee7SLeila Ghaffari   bc->num_wall = 16;
1352fe7aee7SLeila Ghaffari   PetscBool flg;
1362fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_wall",
1372fe7aee7SLeila Ghaffari                               "Face IDs to apply wall BC",
1382fe7aee7SLeila Ghaffari                               NULL, bc->walls, &bc->num_wall, NULL); CHKERRQ(ierr);
1392fe7aee7SLeila Ghaffari   bc->num_comps = 5;
1402fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-wall_comps",
1412fe7aee7SLeila Ghaffari                               "An array of constrained component numbers",
1422fe7aee7SLeila Ghaffari                               NULL, bc->wall_comps, &bc->num_comps, &flg); CHKERRQ(ierr);
1432fe7aee7SLeila Ghaffari   // Slip Boundary Conditions
1442fe7aee7SLeila Ghaffari   for (PetscInt j=0; j<3; j++) {
1452fe7aee7SLeila Ghaffari     bc->num_slip[j] = 16;
1462fe7aee7SLeila Ghaffari     PetscBool flg;
1472fe7aee7SLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1482fe7aee7SLeila Ghaffari     ierr = PetscOptionsIntArray(flags[j],
1492fe7aee7SLeila Ghaffari                                 "Face IDs to apply slip BC",
1502fe7aee7SLeila Ghaffari                                 NULL, bc->slips[j], &bc->num_slip[j], &flg); CHKERRQ(ierr);
1512fe7aee7SLeila Ghaffari     if (flg) bc->user_bc = PETSC_TRUE;
1522fe7aee7SLeila Ghaffari   }
1532fe7aee7SLeila Ghaffari 
1542fe7aee7SLeila Ghaffari   // Error if wall and slip BCs are set on the same face
1552fe7aee7SLeila Ghaffari   if (bc->user_bc)
1562fe7aee7SLeila Ghaffari     for (PetscInt c = 0; c < 3; c++)
1572fe7aee7SLeila Ghaffari       for (PetscInt s = 0; s < bc->num_slip[c]; s++)
1582fe7aee7SLeila Ghaffari         for (PetscInt w = 0; w < bc->num_wall; w++)
1592fe7aee7SLeila Ghaffari           if (bc->slips[c][s] == bc->walls[w])
1607578c821SJeremy L Thompson             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
16108140895SJed Brown                     "Boundary condition already set on face %" PetscInt_FMT "!\n",
1622fe7aee7SLeila Ghaffari                     bc->walls[w]);
1632fe7aee7SLeila Ghaffari 
1642fe7aee7SLeila Ghaffari   // Inflow BCs
1652fe7aee7SLeila Ghaffari   bc->num_inflow = 16;
1662fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_inflow",
1672fe7aee7SLeila Ghaffari                               "Face IDs to apply inflow BC",
1682fe7aee7SLeila Ghaffari                               NULL, bc->inflows, &bc->num_inflow, NULL); CHKERRQ(ierr);
1692fe7aee7SLeila Ghaffari   // Outflow BCs
1702fe7aee7SLeila Ghaffari   bc->num_outflow = 16;
1712fe7aee7SLeila Ghaffari   ierr = PetscOptionsIntArray("-bc_outflow",
1722fe7aee7SLeila Ghaffari                               "Face IDs to apply outflow BC",
1732fe7aee7SLeila Ghaffari                               NULL, bc->outflows, &bc->num_outflow, NULL); CHKERRQ(ierr);
174f4ca79c2SJames Wright   // Freestream BCs
175f4ca79c2SJames Wright   bc->num_freestream = 16;
176f4ca79c2SJames Wright   ierr = PetscOptionsIntArray("-bc_freestream",
177f4ca79c2SJames Wright                               "Face IDs to apply freestream BC",
178f4ca79c2SJames Wright                               NULL, bc->freestreams, &bc->num_freestream, NULL); CHKERRQ(ierr);
1792fe7aee7SLeila Ghaffari 
18067490bc6SJeremy L Thompson   PetscOptionsEnd();
18177841947SLeila Ghaffari 
18277841947SLeila Ghaffari   PetscFunctionReturn(0);
18377841947SLeila Ghaffari }
184