xref: /honee/src/cloptions.c (revision fce2147ef5906fbd4ea41b4f5cee0ce9c3d579ee)
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 
2006f41313SJames Wright   PetscFunctionBeginUser;
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, "blasius", NS_BLASIUS));
262b916ea7SJeremy L Thompson   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "channel", NS_CHANNEL));
27e7754af5SKenneth E. Jansen   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "gaussian_wave", NS_GAUSSIAN_WAVE));
28b8fb7609SAdeleke O. Bankole   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "newtonian", NS_NEWTONIAN_IG));
29692bc0d9SJames Wright   PetscCall(PetscFunctionListAdd(&app_ctx->problems, "taylor_green", NS_TAYLOR_GREEN));
30d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
31a515125bSLeila Ghaffari }
32a515125bSLeila Ghaffari 
33a515125bSLeila Ghaffari // Process general command line options
342b916ea7SJeremy L Thompson PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx, SimpleBC bc) {
35a515125bSLeila Ghaffari   PetscBool ceed_flag    = PETSC_FALSE;
36a515125bSLeila Ghaffari   PetscBool problem_flag = PETSC_FALSE;
3791a36801SJames Wright   PetscBool option_set   = PETSC_FALSE;
382b916ea7SJeremy L Thompson 
39a515125bSLeila Ghaffari   PetscFunctionBeginUser;
402b916ea7SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "Navier-Stokes in PETSc with libCEED", NULL);
41a515125bSLeila Ghaffari 
422b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
432b916ea7SJeremy L Thompson                                sizeof(app_ctx->ceed_resource), &ceed_flag));
44a515125bSLeila Ghaffari 
450e1e9333SJames Wright   app_ctx->test_type = TESTTYPE_NONE;
460e1e9333SJames Wright   PetscCall(PetscOptionsEnum("-test_type", "Type of test to run", NULL, TestTypes, (PetscEnum)(app_ctx->test_type), (PetscEnum *)&app_ctx->test_type,
470e1e9333SJames Wright                              NULL));
48a515125bSLeila Ghaffari 
49a515125bSLeila Ghaffari   app_ctx->test_tol = 1E-11;
502b916ea7SJeremy L Thompson   PetscCall(PetscOptionsScalar("-compare_final_state_atol", "Test absolute tolerance", NULL, app_ctx->test_tol, &app_ctx->test_tol, NULL));
51a515125bSLeila Ghaffari 
52c931fa59SJames Wright   PetscCall(PetscOptionsString("-compare_final_state_filename", "Test filename", NULL, app_ctx->test_file_path, app_ctx->test_file_path,
53c931fa59SJames Wright                                sizeof(app_ctx->test_file_path), NULL));
54a515125bSLeila Ghaffari 
552b916ea7SJeremy L Thompson   PetscCall(PetscOptionsFList("-problem", "Problem to solve", NULL, app_ctx->problems, app_ctx->problem_name, app_ctx->problem_name,
562b916ea7SJeremy L Thompson                               sizeof(app_ctx->problem_name), &problem_flag));
57a515125bSLeila Ghaffari 
58a515125bSLeila Ghaffari   app_ctx->viz_refine = 0;
592b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-viz_refine", "Regular refinement levels for visualization", NULL, app_ctx->viz_refine, &app_ctx->viz_refine, NULL));
60a515125bSLeila Ghaffari 
61852e5969SJed Brown   app_ctx->checkpoint_interval = 10;
62852e5969SJed Brown   app_ctx->checkpoint_vtk      = PETSC_FALSE;
63852e5969SJed Brown   PetscCall(PetscOptionsDeprecated("-output_freq", "-checkpoint_interval", "libCEED 0.11.1", "Use -checkpoint_vtk true to include VTK output"));
64852e5969SJed Brown   PetscCall(PetscOptionsInt("-output_freq", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval,
65852e5969SJed Brown                             &app_ctx->checkpoint_interval, &option_set));
66852e5969SJed Brown   if (option_set) app_ctx->checkpoint_vtk = PETSC_TRUE;
67852e5969SJed Brown   PetscCall(PetscOptionsInt("-checkpoint_interval", "Frequency of output, in number of steps", NULL, app_ctx->checkpoint_interval,
68852e5969SJed Brown                             &app_ctx->checkpoint_interval, NULL));
69852e5969SJed Brown   PetscCall(PetscOptionsBool("-checkpoint_vtk", "Include VTK (*.vtu) output at each binary checkpoint", NULL, app_ctx->checkpoint_vtk,
70852e5969SJed Brown                              &app_ctx->checkpoint_vtk, NULL));
71a515125bSLeila Ghaffari 
722b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-output_add_stepnum2bin", "Add step number to the binary outputs", NULL, app_ctx->add_stepnum2bin,
732b916ea7SJeremy L Thompson                              &app_ctx->add_stepnum2bin, NULL));
7491a36801SJames Wright 
7591a36801SJames Wright   PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2));
762b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL));
7791a36801SJames Wright 
78a515125bSLeila Ghaffari   app_ctx->cont_steps = 0;
792b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-continue", "Continue from previous solution", NULL, app_ctx->cont_steps, &app_ctx->cont_steps, NULL));
80a515125bSLeila Ghaffari 
8191a36801SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_file, "[output_dir]/ns-solution.bin"));
822b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-continue_filename", "Filename to get initial condition from", NULL, app_ctx->cont_file, app_ctx->cont_file,
8391a36801SJames Wright                                sizeof(app_ctx->cont_file), &option_set));
842b916ea7SJeremy L Thompson   if (!option_set) PetscCall(PetscSNPrintf(app_ctx->cont_file, sizeof app_ctx->cont_file, "%s/ns-solution.bin", app_ctx->output_dir));
859293eaa1SJed Brown   if (option_set && app_ctx->cont_steps == 0) app_ctx->cont_steps = -1;  // Read time from file
8691a36801SJames Wright 
8791a36801SJames Wright   PetscCall(PetscStrcpy(app_ctx->cont_time_file, "[output_dir]/ns-time.bin"));
882b916ea7SJeremy L Thompson   PetscCall(PetscOptionsString("-continue_time_filename", "Filename to get initial condition time from", NULL, app_ctx->cont_time_file,
892b916ea7SJeremy L Thompson                                app_ctx->cont_time_file, sizeof(app_ctx->cont_time_file), &option_set));
902b916ea7SJeremy 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));
9191a36801SJames Wright 
92a515125bSLeila Ghaffari   app_ctx->degree = 1;
932b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of finite elements", NULL, app_ctx->degree, &app_ctx->degree, NULL));
94a515125bSLeila Ghaffari 
951219168aSLeila Ghaffari   app_ctx->q_extra = 0;
962b916ea7SJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL));
97a515125bSLeila Ghaffari 
98b107fddaSJed Brown   {
99b107fddaSJed Brown     PetscBool option_set;
100b107fddaSJed Brown     char      amat_type[256] = "";
1012b916ea7SJeremy L Thompson     PetscCall(PetscOptionsFList("-amat_type", "Set the type of Amat distinct from Pmat (-dm_mat_type)", NULL, MatList, amat_type, amat_type,
1022b916ea7SJeremy L Thompson                                 sizeof(amat_type), &option_set));
1032b916ea7SJeremy L Thompson     if (option_set) PetscCall(PetscStrallocpy(amat_type, (char **)&app_ctx->amat_type));
104b107fddaSJed Brown   }
1052b916ea7SJeremy L Thompson   PetscCall(PetscOptionsBool("-pmat_pbdiagonal", "Assemble only point-block diagonal for Pmat", NULL, app_ctx->pmat_pbdiagonal,
106b107fddaSJed Brown                              &app_ctx->pmat_pbdiagonal, NULL));
107b107fddaSJed Brown 
108a515125bSLeila Ghaffari   // Provide default ceed resource if not specified
109a515125bSLeila Ghaffari   if (!ceed_flag) {
110a515125bSLeila Ghaffari     const char *ceed_resource = "/cpu/self";
111a515125bSLeila Ghaffari     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
112a515125bSLeila Ghaffari   }
113a467869fSJed Brown   // If we request a GPU, make sure PETSc has initialized its device (which is
114a467869fSJed Brown   // MPI-aware in case multiple devices are available) before CeedInit so that
115a467869fSJed Brown   // PETSc and libCEED agree about which device to use.
116a467869fSJed Brown   if (strncmp(app_ctx->ceed_resource, "/gpu", 4) == 0) PetscCall(PetscDeviceInitialize(PETSC_DEVICE_DEFAULT()));
117a515125bSLeila Ghaffari 
118a515125bSLeila Ghaffari   // Provide default problem if not specified
119a515125bSLeila Ghaffari   if (!problem_flag) {
120a515125bSLeila Ghaffari     const char *problem_name = "density_current";
121a515125bSLeila Ghaffari     strncpy(app_ctx->problem_name, problem_name, 16);
122a515125bSLeila Ghaffari   }
123a515125bSLeila Ghaffari 
124002797a3SLeila Ghaffari   // Wall Boundary Conditions
125002797a3SLeila Ghaffari   bc->num_wall = 16;
126002797a3SLeila Ghaffari   PetscBool flg;
1272b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_wall", "Face IDs to apply wall BC", NULL, bc->walls, &bc->num_wall, NULL));
128002797a3SLeila Ghaffari   bc->num_comps = 5;
1292b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-wall_comps", "An array of constrained component numbers", NULL, bc->wall_comps, &bc->num_comps, &flg));
130*fce2147eSJames Wright 
131*fce2147eSJames Wright   {  // Symmetry Boundary Conditions
132*fce2147eSJames Wright     const char *deprecated[3] = {"-bc_slip_x", "-bc_slip_y", "-bc_slip_z"};
133*fce2147eSJames Wright     const char *flags[3]      = {"-bc_symmetry_x", "-bc_symmetry_y", "-bc_symmetry_z"};
134*fce2147eSJames Wright     PetscBool   flg, has_symmetry = PETSC_FALSE;
135*fce2147eSJames Wright 
136002797a3SLeila Ghaffari     for (PetscInt j = 0; j < 3; j++) {
137*fce2147eSJames Wright       bc->num_symmetry[j] = 16;
138*fce2147eSJames Wright       PetscCall(PetscOptionsDeprecated(deprecated[j], flags[j], "libCEED 0.12.0",
139*fce2147eSJames Wright                                        "Use -bc_symmetry_[x,y,z] for direct equivalency, or -bc_slip for weak, Riemann-based, direction-invariant "
140*fce2147eSJames Wright                                        "slip/no-penatration boundary conditions"));
141*fce2147eSJames Wright       PetscCall(PetscOptionsIntArray(flags[j], "Face IDs to apply symmetry BC", NULL, bc->symmetries[j], &bc->num_symmetry[j], &flg));
142*fce2147eSJames Wright       if (!flg) {
143*fce2147eSJames Wright         bc->num_symmetry[j] = 16;
144*fce2147eSJames Wright         PetscCall(PetscOptionsIntArray(deprecated[j], "Face IDs to apply slip BC", NULL, bc->symmetries[j], &bc->num_symmetry[j], &flg));
145*fce2147eSJames Wright       }
146*fce2147eSJames Wright       if (bc->num_symmetry[j] > 0) has_symmetry = PETSC_TRUE;
147002797a3SLeila Ghaffari     }
148002797a3SLeila Ghaffari 
149*fce2147eSJames Wright     // Error if wall and symmetry BCs are set on the same face
150*fce2147eSJames Wright     if (has_symmetry) {
151e419654dSJeremy L Thompson       for (PetscInt c = 0; c < 3; c++) {
152*fce2147eSJames Wright         for (PetscInt s = 0; s < bc->num_symmetry[c]; s++) {
153e419654dSJeremy L Thompson           for (PetscInt w = 0; w < bc->num_wall; w++) {
154*fce2147eSJames Wright             PetscCheck(bc->symmetries[c][s] != bc->walls[w], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG,
1555d28dccaSJames Wright                        "Boundary condition already set on face %" PetscInt_FMT "!\n", bc->walls[w]);
156e419654dSJeremy L Thompson           }
157e419654dSJeremy L Thompson         }
158e419654dSJeremy L Thompson       }
159e419654dSJeremy L Thompson     }
160*fce2147eSJames Wright   }
16191933550SJames Wright   app_ctx->wall_forces.num_wall = bc->num_wall;
16291933550SJames Wright   PetscMalloc1(bc->num_wall, &app_ctx->wall_forces.walls);
16391933550SJames Wright   PetscCall(PetscArraycpy(app_ctx->wall_forces.walls, bc->walls, bc->num_wall));
164002797a3SLeila Ghaffari 
165002797a3SLeila Ghaffari   // Inflow BCs
166002797a3SLeila Ghaffari   bc->num_inflow = 16;
1672b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_inflow", "Face IDs to apply inflow BC", NULL, bc->inflows, &bc->num_inflow, NULL));
168002797a3SLeila Ghaffari   // Outflow BCs
169002797a3SLeila Ghaffari   bc->num_outflow = 16;
1702b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_outflow", "Face IDs to apply outflow BC", NULL, bc->outflows, &bc->num_outflow, NULL));
171df55ba5fSJames Wright   // Freestream BCs
172df55ba5fSJames Wright   bc->num_freestream = 16;
1732b916ea7SJeremy L Thompson   PetscCall(PetscOptionsIntArray("-bc_freestream", "Face IDs to apply freestream BC", NULL, bc->freestreams, &bc->num_freestream, NULL));
174002797a3SLeila Ghaffari 
175b0488d1fSJames Wright   // Statistics Options
176c931fa59SJames Wright   app_ctx->turb_spanstats_collect_interval = 1;
177c931fa59SJames Wright   PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_collect_interval", "Number of timesteps between statistics collection", NULL,
178c931fa59SJames Wright                             app_ctx->turb_spanstats_collect_interval, &app_ctx->turb_spanstats_collect_interval, NULL));
179b0488d1fSJames Wright 
180c931fa59SJames Wright   app_ctx->turb_spanstats_viewer_interval = -1;
181c931fa59SJames Wright   PetscCall(PetscOptionsInt("-ts_monitor_turbulence_spanstats_viewer_interval", "Number of timesteps between statistics viewer writing", NULL,
182c931fa59SJames Wright                             app_ctx->turb_spanstats_viewer_interval, &app_ctx->turb_spanstats_viewer_interval, NULL));
183b0488d1fSJames Wright 
184c931fa59SJames Wright   PetscCall(PetscOptionsViewer("-ts_monitor_turbulence_spanstats_viewer", "Viewer for the statistics", NULL, &app_ctx->turb_spanstats_viewer,
185c931fa59SJames Wright                                &app_ctx->turb_spanstats_viewer_format, &app_ctx->turb_spanstats_enable));
186b0488d1fSJames Wright 
187c5e9980aSAdeleke O. Bankole   PetscCall(PetscOptionsViewer("-ts_monitor_wall_force", "Viewer for force on each (no-slip) wall", NULL, &app_ctx->wall_forces.viewer,
188c5e9980aSAdeleke O. Bankole                                &app_ctx->wall_forces.viewer_format, NULL));
189c5e9980aSAdeleke O. Bankole 
19001ab89c1SJames Wright   // SGS Model Options
19101ab89c1SJames Wright   app_ctx->sgs_model_type = SGS_MODEL_NONE;
19201ab89c1SJames Wright   PetscCall(PetscOptionsEnum("-sgs_model_type", "Subgrid Stress Model type", NULL, SGSModelTypes, (PetscEnum)app_ctx->sgs_model_type,
19301ab89c1SJames Wright                              (PetscEnum *)&app_ctx->sgs_model_type, NULL));
19401ab89c1SJames Wright 
19588b07121SJames Wright   PetscCall(PetscOptionsBool("-diff_filter_monitor", "Enable differential filtering TSMonitor", NULL, app_ctx->diff_filter_monitor,
19688b07121SJames Wright                              &app_ctx->diff_filter_monitor, NULL));
19788b07121SJames Wright 
198f31f4833SJames Wright   // Mesh Transformation Options
199f31f4833SJames Wright   app_ctx->mesh_transform_type = MESH_TRANSFORM_NONE;
200f31f4833SJames Wright   PetscCall(PetscOptionsEnum("-mesh_transform", "Mesh transform to perform", NULL, MeshTransformTypes, (PetscEnum)app_ctx->mesh_transform_type,
201f31f4833SJames Wright                              (PetscEnum *)&app_ctx->mesh_transform_type, NULL));
202f31f4833SJames Wright 
2031c17f66aSJames Wright   PetscCall(
2041c17f66aSJames Wright       PetscOptionsBool("-sgs_train_enable", "Enable Data-Driven SGS training", NULL, app_ctx->sgs_train_enable, &app_ctx->sgs_train_enable, NULL));
2051c17f66aSJames Wright 
2061485969bSJeremy L Thompson   PetscOptionsEnd();
207d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
208a515125bSLeila Ghaffari }
209