xref: /honee/src/cloptions.c (revision 5d28dccaccae4dbbdfc8aa7c8439b84e1bbb591b)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Command line option processing for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11a467869fSJed Brown #include <petscdevice.h>
12e419654dSJeremy L Thompson #include <petscsys.h>
13a467869fSJed Brown 
14a515125bSLeila Ghaffari #include "../navierstokes.h"
15a515125bSLeila Ghaffari 
16a515125bSLeila Ghaffari // Register problems to be available on the command line
17a515125bSLeila Ghaffari PetscErrorCode RegisterProblems_NS(AppCtx app_ctx) {
18a515125bSLeila Ghaffari   app_ctx->problems = NULL;
19a515125bSLeila Ghaffari   PetscFunctionBeginUser;
20a515125bSLeila Ghaffari 
212b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "density_current", NS_DENSITY_CURRENT));
222b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "euler_vortex", NS_EULER_VORTEX));
232b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "shocktube", NS_SHOCKTUBE));
242b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection", NS_ADVECTION));
252b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "advection2d", NS_ADVECTION2D));
262b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "blasius", NS_BLASIUS));
272b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "channel", NS_CHANNEL));
28e7754af5SKenneth E. Jansen   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "gaussian_wave", NS_GAUSSIAN_WAVE));
29b8fb7609SAdeleke O. Bankole   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "newtonian", NS_NEWTONIAN_IG));
30bb8a0c61SJames Wright 
31a515125bSLeila Ghaffari   PetscFunctionReturn(0);
32a515125bSLeila Ghaffari }
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari // Process general command line options
352b916ea7SJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc) {
36a515125bSLeila Ghaffari   PetscBool ceed_flag    = PETSC_FALSE;
37a515125bSLeila Ghaffari   PetscBool problem_flag = PETSC_FALSE;
3891a36801SJames Wright   PetscBool option_set   = PETSC_FALSE;
392b916ea7SJeremy L Thompson 
40a515125bSLeila Ghaffari   PetscFunctionBeginUser;
41a515125bSLeila Ghaffari 
422b916ea7SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", NULL);
43a515125bSLeila Ghaffari 
442b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
452b916ea7SJeremy L Thompson                                sizeof(app_ctx->ceed_resource), &ceed_flag));
46a515125bSLeila Ghaffari 
470e1e9333SJames Wright   app_ctx->test_type = TESTTYPE_NONE;
480e1e9333SJames Wright   PetscCall(PetscOptionsEnum("-test_type", "Type of test to run", NULL, TestTypes, (PetscEnum)(app_ctx->test_type), (PetscEnum *)&app_ctx->test_type,
490e1e9333SJames Wright                              NULL));
50a515125bSLeila Ghaffari 
51a515125bSLeila Ghaffari   app_ctx->test_tol = 1E-11;
522b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-compare_final_state_atol", "Test absolute tolerance", NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL));
53a515125bSLeila Ghaffari 
54c931fa59SJames Wright   PetscCall(PetscOptionsString("-compare_final_state_filename", "Test filename", NULL, app_ctx->test_file_path, app_ctx->test_file_path,
55c931fa59SJames Wright                                sizeof(app_ctx->test_file_path), NULL));
56a515125bSLeila Ghaffari 
572b916ea7SJeremy L Thompson   PetscCall(PetscOptionsFList("-problem", "Problem to solve", NULL, app_ctx->problems, app_ctx->problem_name, app_ctx->problem_name,
582b916ea7SJeremy L Thompson                               sizeof(app_ctx->problem_name), &problem_flag));
59a515125bSLeila Ghaffari 
60a515125bSLeila Ghaffari   app_ctx->viz_refine = 0;
612b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-viz_refine", "Regular refinement levels for visualization", NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL));
62a515125bSLeila Ghaffari 
63852e5969SJed Brown   app_ctx->checkpoint_interval = 10;
64852e5969SJed Brown   app_ctx->checkpoint_vtk      = PETSC_FALSE;
65852e5969SJed Brown   PetscCall(PetscOptionsDeprecated("-output_freq", "-checkpoint_interval", "libCEED 0.11.1", "Use -checkpoint_vtk true to include VTK output"));
66852e5969SJed Brown   PetscCall(PetscOptionsInt("-output_freq", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval,
67852e5969SJed Brown                             &app_ctx->checkpoint_interval, &option_set));
68852e5969SJed Brown   if (option_set) app_ctx->checkpoint_vtk = PETSC_TRUE;
69852e5969SJed Brown   PetscCall(PetscOptionsInt("-checkpoint_interval", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval,
70852e5969SJed Brown                             &app_ctx->checkpoint_interval, NULL));
71852e5969SJed Brown   PetscCall(PetscOptionsBool("-checkpoint_vtk", "Include VTK (*.vtu) output at each binary checkpoint", NULL, app_ctx->checkpoint_vtk,
72852e5969SJed Brown                              &app_ctx->checkpoint_vtk, NULL));
73a515125bSLeila Ghaffari 
742b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-output_add_stepnum2bin", "Add step number to the binary outputs", NULL, app_ctx->add_stepnum2bin,
752b916ea7SJeremy L Thompson                              &app_ctx->add_stepnum2bin, NULL));
7691a36801SJames Wright 
7791a36801SJames Wright   PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2));
782b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL));
7991a36801SJames Wright 
80a515125bSLeila Ghaffari   app_ctx->cont_steps = 0;
812b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-continue", "Continue from previous solution", NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL));
82a515125bSLeila Ghaffari 
8391a36801SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_file, "[output_dir]/ns-solution.bin"));
842b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-continue_filename", "Filename to get initial condition from", NULL, app_ctx->cont_file, app_ctx->cont_file,
8591a36801SJames Wright                                sizeof(app_ctx->cont_file), &option_set));
862b916ea7SJeremy L Thompson   if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_file, sizeof app_ctx->cont_file, "%s/ns-solution.bin", app_ctx->output_dir));
879293eaa1SJed Brown   if (option_set && app_ctx->cont_steps == 0) app_ctx->cont_steps = -1;  // Read time from file
8891a36801SJames Wright 
8991a36801SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_time_file, "[output_dir]/ns-time.bin"));
902b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-continue_time_filename", "Filename to get initial condition time from", NULL, app_ctx->cont_time_file,
912b916ea7SJeremy L Thompson                                app_ctx->cont_time_file, sizeof(app_ctx->cont_time_file), &option_set));
922b916ea7SJeremy 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));
9391a36801SJames Wright 
94a515125bSLeila Ghaffari   app_ctx->degree = 1;
952b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of finite elements", NULL, app_ctx->degree, &app_ctx->degree, NULL));
96a515125bSLeila Ghaffari 
971219168aSLeila Ghaffari   app_ctx->q_extra = 0;
982b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL));
99a515125bSLeila Ghaffari 
100b107fddaSJed Brown   {
101b107fddaSJed Brown     PetscBool option_set;
102b107fddaSJed Brown     char      amat_type[256] = "";
1032b916ea7SJeremy L Thompson     PetscCall(PetscOptionsFList("-amat_type", "Set the type of Amat distinct from Pmat (-dm_mat_type)", NULL, MatList, amat_type, amat_type,
1042b916ea7SJeremy L Thompson                                 sizeof(amat_type), &option_set));
1052b916ea7SJeremy L Thompson     if (option_set) PetscCall(PetscStrallocpy(amat_type, (char **)&app_ctx->amat_type));
106b107fddaSJed Brown   }
1072b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-pmat_pbdiagonal", "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal,
108b107fddaSJed Brown                              &app_ctx->pmat_pbdiagonal, NULL));
109b107fddaSJed Brown 
110a515125bSLeila Ghaffari   // Provide default ceed resource if not specified
111a515125bSLeila Ghaffari   if (!ceed_flag) {
112a515125bSLeila Ghaffari     const char *ceed_resource = "/cpu/self";
113a515125bSLeila Ghaffari     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
114a515125bSLeila Ghaffari   }
115a467869fSJed Brown   // If we request a GPU, make sure PETSc has initialized its device (which is
116a467869fSJed Brown   // MPI-aware in case multiple devices are available) before CeedInit so that
117a467869fSJed Brown   // PETSc and libCEED agree about which device to use.
118a467869fSJed Brown   if (strncmp(app_ctx->ceed_resource, "/gpu", 4) == 0) PetscCall(PetscDeviceInitialize(PETSC_DEVICE_DEFAULT()));
119a515125bSLeila Ghaffari 
120a515125bSLeila Ghaffari   // Provide default problem if not specified
121a515125bSLeila Ghaffari   if (!problem_flag) {
122a515125bSLeila Ghaffari     const char *problem_name = "density_current";
123a515125bSLeila Ghaffari     strncpy(app_ctx->problem_name, problem_name, 16);
124a515125bSLeila Ghaffari   }
125a515125bSLeila Ghaffari 
126002797a3SLeila Ghaffari   // Wall Boundary Conditions
127002797a3SLeila Ghaffari   bc->num_wall = 16;
128002797a3SLeila Ghaffari   PetscBool flg;
1292b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_wall", "Face IDs to apply wall BC", NULL, bc->walls, &bc->num_wall, NULL));
130002797a3SLeila Ghaffari   bc->num_comps = 5;
1312b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-wall_comps", "An array of constrained component numbers", NULL, bc->wall_comps, &bc->num_comps, &flg));
132002797a3SLeila Ghaffari   // Slip Boundary Conditions
133002797a3SLeila Ghaffari   for (PetscInt j = 0; j < 3; j++) {
134002797a3SLeila Ghaffari     bc->num_slip[j] = 16;
135002797a3SLeila Ghaffari     PetscBool   flg;
136002797a3SLeila Ghaffari     const char *flags[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
1372b916ea7SJeremy L Thompson     PetscCall(PetscOptionsIntArray(flags[j], "Face IDs to apply slip BC", NULL, bc->slips[j], &bc->num_slip[j], &flg));
138002797a3SLeila Ghaffari     if (flg) bc->user_bc = PETSC_TRUE;
139002797a3SLeila Ghaffari   }
140002797a3SLeila Ghaffari 
141002797a3SLeila Ghaffari   // Error if wall and slip BCs are set on the same face
142e419654dSJeremy L Thompson   if (bc->user_bc) {
143e419654dSJeremy L Thompson     for (PetscInt c = 0; c < 3; c++) {
144e419654dSJeremy L Thompson       for (PetscInt s = 0; s < bc->num_slip[c]; s++) {
145e419654dSJeremy L Thompson         for (PetscInt w = 0; w < bc->num_wall; w++) {
146*5d28dccaSJames Wright           PetscCheck(bc->slips[c][s] != bc->walls[w], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
147*5d28dccaSJames Wright                      "Boundary condition already set on face %" PetscInt_FMT "!\n", bc->walls[w]);
148e419654dSJeremy L Thompson         }
149e419654dSJeremy L Thompson       }
150e419654dSJeremy L Thompson     }
151e419654dSJeremy L Thompson   }
152002797a3SLeila Ghaffari 
153002797a3SLeila Ghaffari   // Inflow BCs
154002797a3SLeila Ghaffari   bc->num_inflow = 16;
1552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_inflow", "Face IDs to apply inflow BC", NULL, bc->inflows, &bc->num_inflow, NULL));
156002797a3SLeila Ghaffari   // Outflow BCs
157002797a3SLeila Ghaffari   bc->num_outflow = 16;
1582b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_outflow", "Face IDs to apply outflow BC", NULL, bc->outflows, &bc->num_outflow, NULL));
159df55ba5fSJames Wright   // Freestream BCs
160df55ba5fSJames Wright   bc->num_freestream = 16;
1612b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_freestream", "Face IDs to apply freestream BC", NULL, bc->freestreams, &bc->num_freestream, NULL));
162002797a3SLeila Ghaffari 
163b0488d1fSJames Wright   // Statistics Options
164c931fa59SJames Wright   app_ctx->turb_spanstats_collect_interval = 1;
165c931fa59SJames Wright   PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_collect_interval", "Number of timesteps between statistics collection", NULL,
166c931fa59SJames Wright                             app_ctx->turb_spanstats_collect_interval, &app_ctx->turb_spanstats_collect_interval, NULL));
167b0488d1fSJames Wright 
168c931fa59SJames Wright   app_ctx->turb_spanstats_viewer_interval = -1;
169c931fa59SJames Wright   PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_viewer_interval", "Number of timesteps between statistics viewer writing", NULL,
170c931fa59SJames Wright                             app_ctx->turb_spanstats_viewer_interval, &app_ctx->turb_spanstats_viewer_interval, NULL));
171b0488d1fSJames Wright 
172c931fa59SJames Wright   PetscCall(PetscOptionsViewer("-ts_monitor_turbulence_spanstats_viewer", "Viewer for the statistics", NULL, &app_ctx->turb_spanstats_viewer,
173c931fa59SJames Wright                                &app_ctx->turb_spanstats_viewer_format, &app_ctx->turb_spanstats_enable));
174b0488d1fSJames Wright 
175c5e9980aSAdeleke O. Bankole   PetscCall(PetscOptionsViewer("-ts_monitor_wall_force", "Viewer for force on each (no-slip) wall", NULL, &app_ctx->wall_forces.viewer,
176c5e9980aSAdeleke O. Bankole                                &app_ctx->wall_forces.viewer_format, NULL));
177c5e9980aSAdeleke O. Bankole 
17801ab89c1SJames Wright   // SGS Model Options
17901ab89c1SJames Wright   app_ctx->sgs_model_type = SGS_MODEL_NONE;
18001ab89c1SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_type", "Subgrid Stress Model type", NULL, SGSModelTypes, (PetscEnum)app_ctx->sgs_model_type,
18101ab89c1SJames Wright                              (PetscEnum *)&app_ctx->sgs_model_type, NULL));
18201ab89c1SJames Wright 
1831485969bSJeremy L Thompson   PetscOptionsEnd();
184a515125bSLeila Ghaffari 
185a515125bSLeila Ghaffari   PetscFunctionReturn(0);
186a515125bSLeila Ghaffari }
187