xref: /libCEED/examples/solids/src/cl-options.c (revision 2eb0be0b68ebb7cb25cc038250732c3239325ad2)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 /// @file
9 /// Command line option processing for solid mechanics example using PETSc
10 
11 #include "../include/cl-options.h"
12 
13 // -----------------------------------------------------------------------------
14 // Process command line options
15 // -----------------------------------------------------------------------------
16 // Process general command line options
17 PetscErrorCode ProcessCommandLineOptions(MPI_Comm comm, AppCtx app_ctx) {
18   PetscBool ceed_flag = PETSC_FALSE;
19 
20   PetscFunctionBeginUser;
21 
22   PetscOptionsBegin(comm, NULL, "Elasticity / Hyperelasticity in PETSc with libCEED", NULL);
23 
24   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, app_ctx->ceed_resource, app_ctx->ceed_resource,
25                                sizeof(app_ctx->ceed_resource), &ceed_flag));
26 
27   PetscCall(PetscStrncpy(app_ctx->output_dir, ".", 2));  // Default - current directory
28   PetscCall(PetscOptionsString("-output_dir", "Output directory", NULL, app_ctx->output_dir, app_ctx->output_dir, sizeof(app_ctx->output_dir), NULL));
29 
30   app_ctx->degree = 3;
31   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, app_ctx->degree, &app_ctx->degree, NULL));
32 
33   app_ctx->q_extra = 0;
34   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, app_ctx->q_extra, &app_ctx->q_extra, NULL));
35 
36   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, app_ctx->mesh_file, app_ctx->mesh_file, sizeof(app_ctx->mesh_file), NULL));
37 
38   app_ctx->problem_choice = ELAS_LINEAR;  // Default - Linear Elasticity
39   PetscCall(PetscOptionsEnum("-problem", "Solves Elasticity & Hyperelasticity Problems", NULL, problemTypes, (PetscEnum)app_ctx->problem_choice,
40                              (PetscEnum *)&app_ctx->problem_choice, NULL));
41   app_ctx->name          = problemTypes[app_ctx->problem_choice];
42   app_ctx->name_for_disp = problemTypesForDisp[app_ctx->problem_choice];
43 
44   app_ctx->num_increments = app_ctx->problem_choice == ELAS_LINEAR ? 1 : 10;
45   PetscCall(PetscOptionsInt("-num_steps", "Number of pseudo-time steps", NULL, app_ctx->num_increments, &app_ctx->num_increments, NULL));
46 
47   app_ctx->forcing_choice = FORCE_NONE;  // Default - no forcing term
48   PetscCall(PetscOptionsEnum("-forcing", "Set forcing function option", NULL, forcing_types, (PetscEnum)app_ctx->forcing_choice,
49                              (PetscEnum *)&app_ctx->forcing_choice, NULL));
50 
51   PetscInt max_n             = 3;
52   app_ctx->forcing_vector[0] = 0;
53   app_ctx->forcing_vector[1] = -1;
54   app_ctx->forcing_vector[2] = 0;
55   PetscCall(PetscOptionsScalarArray("-forcing_vec", "Direction to apply constant force", NULL, app_ctx->forcing_vector, &max_n, NULL));
56 
57   if ((app_ctx->problem_choice == ELAS_FSInitial_NH1 || app_ctx->problem_choice == ELAS_FSInitial_NH2 ||
58        app_ctx->problem_choice == ELAS_FSCurrent_NH1 || app_ctx->problem_choice == ELAS_FSCurrent_NH2 ||
59        app_ctx->problem_choice == ELAS_FSInitial_MR1) &&
60       app_ctx->forcing_choice == FORCE_CONST)
61     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,
62             "Cannot use constant forcing and finite strain formulation. "
63             "Constant forcing in reference frame currently unavaliable.");
64 
65   // Dirichlet boundary conditions
66   app_ctx->bc_clamp_count = 16;
67   PetscCall(
68       PetscOptionsIntArray("-bc_clamp", "Face IDs to apply incremental Dirichlet BC", NULL, app_ctx->bc_clamp_faces, &app_ctx->bc_clamp_count, NULL));
69   // Set vector for each clamped BC
70   for (PetscInt i = 0; i < app_ctx->bc_clamp_count; i++) {
71     // Translation vector
72     char         option_name[25];
73     const size_t nclamp_params = sizeof(app_ctx->bc_clamp_max[0]) / sizeof(app_ctx->bc_clamp_max[0][0]);
74     for (PetscInt j = 0; j < nclamp_params; j++) app_ctx->bc_clamp_max[i][j] = 0.;
75 
76     snprintf(option_name, sizeof option_name, "-bc_clamp_%" PetscInt_FMT "_translate", app_ctx->bc_clamp_faces[i]);
77     max_n = 3;
78     PetscCall(PetscOptionsScalarArray(option_name, "Vector to translate clamped end by", NULL, app_ctx->bc_clamp_max[i], &max_n, NULL));
79 
80     // Rotation vector
81     max_n = 5;
82     snprintf(option_name, sizeof option_name, "-bc_clamp_%" PetscInt_FMT "_rotate", app_ctx->bc_clamp_faces[i]);
83     PetscCall(PetscOptionsScalarArray(option_name, "Vector with axis of rotation and rotation, in radians", NULL, &app_ctx->bc_clamp_max[i][3],
84                                       &max_n, NULL));
85 
86     // Normalize
87     PetscScalar norm = sqrt(app_ctx->bc_clamp_max[i][3] * app_ctx->bc_clamp_max[i][3] + app_ctx->bc_clamp_max[i][4] * app_ctx->bc_clamp_max[i][4] +
88                             app_ctx->bc_clamp_max[i][5] * app_ctx->bc_clamp_max[i][5]);
89     if (fabs(norm) < 1e-16) norm = 1;
90     for (PetscInt j = 0; j < 3; j++) app_ctx->bc_clamp_max[i][3 + j] /= norm;
91   }
92 
93   // Neumann boundary conditions
94   app_ctx->bc_traction_count = 16;
95   PetscCall(PetscOptionsIntArray("-bc_traction", "Face IDs to apply traction (Neumann) BC", NULL, app_ctx->bc_traction_faces,
96                                  &app_ctx->bc_traction_count, NULL));
97   // Set vector for each traction BC
98   for (PetscInt i = 0; i < app_ctx->bc_traction_count; i++) {
99     // Translation vector
100     char option_name[25];
101     for (PetscInt j = 0; j < 3; j++) app_ctx->bc_traction_vector[i][j] = 0.;
102 
103     snprintf(option_name, sizeof option_name, "-bc_traction_%" PetscInt_FMT, app_ctx->bc_traction_faces[i]);
104     max_n         = 3;
105     PetscBool set = false;
106     PetscCall(PetscOptionsScalarArray(option_name, "Traction vector for constrained face", NULL, app_ctx->bc_traction_vector[i], &max_n, &set));
107 
108     if (!set) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Traction vector must be set for all traction boundary conditions.");
109   }
110 
111   app_ctx->multigrid_choice = MULTIGRID_LOGARITHMIC;
112   PetscCall(PetscOptionsEnum("-multigrid", "Set multigrid type option", NULL, multigrid_types, (PetscEnum)app_ctx->multigrid_choice,
113                              (PetscEnum *)&app_ctx->multigrid_choice, NULL));
114 
115   app_ctx->test_mode = PETSC_FALSE;
116   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, app_ctx->test_mode, &(app_ctx->test_mode), NULL));
117 
118   app_ctx->expect_final_strain = -1.;
119   PetscCall(PetscOptionsReal("-expect_final_strain_energy", "Expect final strain energy close to this value.", NULL, app_ctx->expect_final_strain,
120                              &app_ctx->expect_final_strain, NULL));
121 
122   app_ctx->test_tol = 1e-8;
123   PetscCall(PetscOptionsReal("-expect_final_state_rtol", "Relative tolerance for final strain energy test", NULL, app_ctx->test_tol,
124                              &app_ctx->test_tol, NULL));
125 
126   app_ctx->view_soln = PETSC_FALSE;
127   PetscCall(PetscOptionsBool("-view_soln", "Write out solution vector for viewing", NULL, app_ctx->view_soln, &(app_ctx->view_soln), NULL));
128 
129   app_ctx->view_final_soln = PETSC_FALSE;
130   PetscCall(PetscOptionsBool("-view_final_soln", "Write out final solution vector for viewing", NULL, app_ctx->view_final_soln,
131                              &(app_ctx->view_final_soln), NULL));
132 
133   PetscBool set;
134   char      energy_viewer_filename[PETSC_MAX_PATH_LEN] = "";
135   PetscCall(PetscOptionsString("-strain_energy_monitor", "Print out current strain energy at every load increment", NULL, energy_viewer_filename,
136                                energy_viewer_filename, sizeof(energy_viewer_filename), &set));
137   if (set) {
138     PetscCall(PetscViewerASCIIOpen(comm, energy_viewer_filename, &app_ctx->energy_viewer));
139     PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "increment,energy\n"));
140     // Initial configuration is base energy state; this may not be true if we extend in the future to initially loaded configurations (because a truly
141     // at-rest initial state may not be realizable).
142     PetscCall(PetscViewerASCIIPrintf(app_ctx->energy_viewer, "%f,%e\n", 0., 0.));
143   }
144   PetscOptionsEnd();  // End of setting AppCtx
145 
146   // Check for all required values set
147   if (app_ctx->test_mode) {
148     if (app_ctx->forcing_choice == FORCE_NONE && !app_ctx->bc_clamp_count) app_ctx->forcing_choice = FORCE_MMS;
149   }
150   if (!app_ctx->bc_clamp_count && app_ctx->forcing_choice != FORCE_MMS) {
151     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "-boundary options needed");
152   }
153 
154   // Provide default ceed resource if not specified
155   if (!ceed_flag) {
156     const char *ceed_resource = "/cpu/self";
157     strncpy(app_ctx->ceed_resource, ceed_resource, 10);
158   }
159 
160   // Determine number of levels
161   switch (app_ctx->multigrid_choice) {
162     case MULTIGRID_LOGARITHMIC:
163       app_ctx->num_levels = ceil(log(app_ctx->degree) / log(2)) + 1;
164       break;
165     case MULTIGRID_UNIFORM:
166       app_ctx->num_levels = app_ctx->degree;
167       break;
168     case MULTIGRID_NONE:
169       app_ctx->num_levels = 1;
170       break;
171   }
172 
173   // Populate array of degrees for each level for multigrid
174   PetscCall(PetscMalloc1(app_ctx->num_levels, &(app_ctx->level_degrees)));
175 
176   switch (app_ctx->multigrid_choice) {
177     case MULTIGRID_LOGARITHMIC:
178       for (int i = 0; i < app_ctx->num_levels - 1; i++) app_ctx->level_degrees[i] = pow(2, i);
179       app_ctx->level_degrees[app_ctx->num_levels - 1] = app_ctx->degree;
180       break;
181     case MULTIGRID_UNIFORM:
182       for (int i = 0; i < app_ctx->num_levels; i++) app_ctx->level_degrees[i] = i + 1;
183       break;
184     case MULTIGRID_NONE:
185       app_ctx->level_degrees[0] = app_ctx->degree;
186       break;
187   }
188 
189   PetscFunctionReturn(0);
190 };
191