xref: /libCEED/examples/fluids/src/cloptions.c (revision 0992c3d11e5379bfb326e16b4ff4fcef892605e9)
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 
11*0992c3d1SJed Brown #include <petscdevice.h>
12*0992c3d1SJed Brown 
1377841947SLeila Ghaffari #include "../navierstokes.h"
1477841947SLeila Ghaffari 
1577841947SLeila Ghaffari // Register problems to be available on the command line
1677841947SLeila Ghaffari PetscErrorCode RegisterProblems_NS(AppCtx app_ctx) {
1777841947SLeila Ghaffari   app_ctx->problems = NULL;
1877841947SLeila Ghaffari   PetscFunctionBeginUser;
1977841947SLeila Ghaffari 
202b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "density_current", NS_DENSITY_CURRENT));
212b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "euler_vortex", NS_EULER_VORTEX));
222b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "shocktube", NS_SHOCKTUBE));
232b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection", NS_ADVECTION));
242b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection2d", NS_ADVECTION2D));
252b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "blasius", NS_BLASIUS));
262b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "channel", NS_CHANNEL));
272b730f8bSJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "newtonian_wave", NS_NEWTONIAN_WAVE));
2888626eedSJames Wright 
2977841947SLeila Ghaffari   PetscFunctionReturn(0);
3077841947SLeila Ghaffari }
3177841947SLeila Ghaffari 
3277841947SLeila Ghaffari // Process general command line options
332b730f8bSJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc) {
3477841947SLeila Ghaffari   PetscBool ceed_flag    = PETSC_FALSE;
3577841947SLeila Ghaffari   PetscBool problem_flag = PETSC_FALSE;
3669293791SJames Wright   PetscBool option_set   = PETSC_FALSE;
372b730f8bSJeremy L Thompson 
3877841947SLeila Ghaffari   PetscFunctionBeginUser;
3977841947SLeila Ghaffari 
402b730f8bSJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", NULL);
4177841947SLeila Ghaffari 
422b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
432b730f8bSJeremy L Thompson                                sizeof(app_ctx->ceed_resource), &ceed_flag));
4477841947SLeila Ghaffari 
4577841947SLeila Ghaffari   app_ctx->test_mode = PETSC_FALSE;
462b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Run in test mode", NULL, app_ctx->test_mode, &app_ctx->test_mode, NULL));
4777841947SLeila Ghaffari 
4877841947SLeila Ghaffari   app_ctx->test_tol = 1E-11;
492b730f8bSJeremy L Thompson   PetscCall(PetscOptionsScalar("-compare_final_state_atol", "Test absolute tolerance", NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL));
5077841947SLeila Ghaffari 
512b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-compare_final_state_filename", "Test filename", NULL, app_ctx->file_path, app_ctx->file_path,
522b730f8bSJeremy L Thompson                                sizeof(app_ctx->file_path), NULL));
5377841947SLeila Ghaffari 
542b730f8bSJeremy L Thompson   PetscCall(PetscOptionsFList("-problem", "Problem to solve", NULL, app_ctx->problems, app_ctx->problem_name, app_ctx->problem_name,
552b730f8bSJeremy L Thompson                               sizeof(app_ctx->problem_name), &problem_flag));
5677841947SLeila Ghaffari 
5777841947SLeila Ghaffari   app_ctx->viz_refine = 0;
582b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-viz_refine", "Regular refinement levels for visualization", NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL));
5977841947SLeila Ghaffari 
6077841947SLeila Ghaffari   app_ctx->output_freq = 10;
612b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-output_freq", "Frequency of output, in number of steps", NULL, app_ctx->output_freq, &app_ctx->output_freq, NULL));
6277841947SLeila Ghaffari 
632b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-output_add_stepnum2bin", "Add step number to the binary outputs", NULL, app_ctx->add_stepnum2bin,
642b730f8bSJeremy L Thompson                              &app_ctx->add_stepnum2bin, NULL));
6569293791SJames Wright 
6669293791SJames Wright   PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2));
672b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL));
6869293791SJames Wright 
6977841947SLeila Ghaffari   app_ctx->cont_steps = 0;
702b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-continue", "Continue from previous solution", NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL));
7177841947SLeila Ghaffari 
7269293791SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_file, "[output_dir]/ns-solution.bin"));
732b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-continue_filename", "Filename to get initial condition from", NULL, app_ctx->cont_file, app_ctx->cont_file,
7469293791SJames Wright                                sizeof(app_ctx->cont_file), &option_set));
752b730f8bSJeremy L Thompson   if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_file, sizeof app_ctx->cont_file, "%s/ns-solution.bin", app_ctx->output_dir));
7669293791SJames Wright 
7769293791SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_time_file, "[output_dir]/ns-time.bin"));
782b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-continue_time_filename", "Filename to get initial condition time from", NULL, app_ctx->cont_time_file,
792b730f8bSJeremy L Thompson                                app_ctx->cont_time_file, sizeof(app_ctx->cont_time_file), &option_set));
802b730f8bSJeremy 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));
8169293791SJames Wright 
8277841947SLeila Ghaffari   app_ctx->degree = 1;
832b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of finite elements", NULL, app_ctx->degree, &app_ctx->degree, NULL));
8477841947SLeila Ghaffari 
8577841947SLeila Ghaffari   app_ctx->q_extra = 2;
862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL));
8777841947SLeila Ghaffari 
88544be873SJed Brown   {
89544be873SJed Brown     PetscBool option_set;
90544be873SJed Brown     char      amat_type[256] = "";
912b730f8bSJeremy L Thompson     PetscCall(PetscOptionsFList("-amat_type", "Set the type of Amat distinct from Pmat (-dm_mat_type)", NULL, MatList, amat_type, amat_type,
922b730f8bSJeremy L Thompson                                 sizeof(amat_type), &option_set));
932b730f8bSJeremy L Thompson     if (option_set) PetscCall(PetscStrallocpy(amat_type, (char **)&app_ctx->amat_type));
94544be873SJed Brown   }
952b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-pmat_pbdiagonal", "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal,
96544be873SJed Brown                              &app_ctx->pmat_pbdiagonal, NULL));
97544be873SJed Brown 
9877841947SLeila Ghaffari   // Provide default ceed resource if not specified
9977841947SLeila Ghaffari   if (!ceed_flag) {
10077841947SLeila Ghaffari     const char *ceed_resource = "/cpu/self";
10177841947SLeila Ghaffari     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
10277841947SLeila Ghaffari   }
103*0992c3d1SJed Brown   // If we request a GPU, make sure PETSc has initialized its device (which is
104*0992c3d1SJed Brown   // MPI-aware in case multiple devices are available) before CeedInit so that
105*0992c3d1SJed Brown   // PETSc and libCEED agree about which device to use.
106*0992c3d1SJed Brown   if (strncmp(app_ctx->ceed_resource, "/gpu", 4) == 0) PetscCall(PetscDeviceInitialize(PETSC_DEVICE_DEFAULT()));
10777841947SLeila Ghaffari 
10877841947SLeila Ghaffari   // Provide default problem if not specified
10977841947SLeila Ghaffari   if (!problem_flag) {
11077841947SLeila Ghaffari     const char *problem_name = "density_current";
11177841947SLeila Ghaffari     strncpy(app_ctx->problem_name, problem_name, 16);
11277841947SLeila Ghaffari   }
11377841947SLeila Ghaffari 
1142fe7aee7SLeila Ghaffari   // Wall Boundary Conditions
1152fe7aee7SLeila Ghaffari   bc->num_wall = 16;
1162fe7aee7SLeila Ghaffari   PetscBool flg;
1172b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_wall", "Face IDs to apply wall BC", NULL, bc->walls, &bc->num_wall, NULL));
1182fe7aee7SLeila Ghaffari   bc->num_comps = 5;
1192b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-wall_comps", "An array of constrained component numbers", NULL, bc->wall_comps, &bc->num_comps, &flg));
1202fe7aee7SLeila Ghaffari   // Slip Boundary Conditions
1212fe7aee7SLeila Ghaffari   for (PetscInt j = 0; j < 3; j++) {
1222fe7aee7SLeila Ghaffari     bc->num_slip[j] = 16;
1232fe7aee7SLeila Ghaffari     PetscBool   flg;
1242fe7aee7SLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1252b730f8bSJeremy L Thompson     PetscCall(PetscOptionsIntArray(flags[j], "Face IDs to apply slip BC", NULL, bc->slips[j], &bc->num_slip[j], &flg));
1262fe7aee7SLeila Ghaffari     if (flg) bc->user_bc = PETSC_TRUE;
1272fe7aee7SLeila Ghaffari   }
1282fe7aee7SLeila Ghaffari 
1292fe7aee7SLeila Ghaffari   // Error if wall and slip BCs are set on the same face
1302fe7aee7SLeila Ghaffari   if (bc->user_bc)
1312fe7aee7SLeila Ghaffari     for (PetscInt c = 0; c < 3; c++)
1322fe7aee7SLeila Ghaffari       for (PetscInt s = 0; s < bc->num_slip[c]; s++)
1332fe7aee7SLeila Ghaffari         for (PetscInt w = 0; w < bc->num_wall; w++)
1342fe7aee7SLeila Ghaffari           if (bc->slips[c][s] == bc->walls[w])
1352b730f8bSJeremy L Thompson             SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Boundary condition already set on face %" PetscInt_FMT "!\n", bc->walls[w]);
1362fe7aee7SLeila Ghaffari 
1372fe7aee7SLeila Ghaffari   // Inflow BCs
1382fe7aee7SLeila Ghaffari   bc->num_inflow = 16;
1392b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_inflow", "Face IDs to apply inflow BC", NULL, bc->inflows, &bc->num_inflow, NULL));
1402fe7aee7SLeila Ghaffari   // Outflow BCs
1412fe7aee7SLeila Ghaffari   bc->num_outflow = 16;
1422b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_outflow", "Face IDs to apply outflow BC", NULL, bc->outflows, &bc->num_outflow, NULL));
143f4ca79c2SJames Wright   // Freestream BCs
144f4ca79c2SJames Wright   bc->num_freestream = 16;
1452b730f8bSJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_freestream", "Face IDs to apply freestream BC", NULL, bc->freestreams, &bc->num_freestream, NULL));
1462fe7aee7SLeila Ghaffari 
14767490bc6SJeremy L Thompson   PetscOptionsEnd();
14877841947SLeila Ghaffari 
14977841947SLeila Ghaffari   PetscFunctionReturn(0);
15077841947SLeila Ghaffari }
151