-
Notifications
You must be signed in to change notification settings - Fork 0
/
parameter-namespace.cc
221 lines (192 loc) · 7.75 KB
/
parameter-namespace.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
#include <deal.II/base/parameter_handler.h>
#include "parameter-namespace.h"
namespace Parameters
{
void FESystem::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Finite element system");
{
prm.declare_entry("Polynomial degree", "2",
Patterns::Integer(0),
"Displacement system polynomial order");
prm.declare_entry("Quadrature order", "3",
Patterns::Integer(0),
"Gauss quadrature order");
}
prm.leave_subsection();
}
void FESystem::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Finite element system");
{
poly_degree = prm.get_integer("Polynomial degree");
quad_order = prm.get_integer("Quadrature order");
}
prm.leave_subsection();
}
// @sect4{Geometry} =======================================================================================
void Geometry::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Geometry");
{
prm.declare_entry("Global refinement", "2",
Patterns::Integer(0),
"Global refinement level");
prm.declare_entry("Grid scale", "1e-3",
Patterns::Double(0.0),
"Global grid scaling factor");
prm.declare_entry("Pressure ratio p/p0", "100",
Patterns::Selection("20|40|60|80|100"),
"Ratio of applied pressure to reference pressure");
}
prm.leave_subsection();
}
void Geometry::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Geometry");
{
global_refinement = prm.get_integer("Global refinement");
scale = prm.get_double("Grid scale");
p_p0 = prm.get_double("Pressure ratio p/p0");
}
prm.leave_subsection();
}
// @sect4{Materials} =============================================================================
void Materials::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Material properties");
{
prm.declare_entry("Poisson's ratio", "0.4999",
Patterns::Double(-1.0,0.5),
"Poisson's ratio");
prm.declare_entry("Shear modulus", "80.194e6",
Patterns::Double(),
"Shear modulus");
}
prm.leave_subsection();
}
void Materials::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Material properties");
{
nu = prm.get_double("Poisson's ratio");
mu = prm.get_double("Shear modulus");
}
prm.leave_subsection();
}
// @sect4{Linear solver} ================================================================================
void LinearSolver::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Linear solver");
{
prm.declare_entry("Solver type", "cg",
Patterns::Selection("cg|Direct"),
"Type of solver used to solve the linear system");
prm.declare_entry("Residual", "1e-6",
Patterns::Double(0.0),
"Linear solver residual (scaled by residual norm)");
prm.declare_entry("Max iteration multiplier", "1",
Patterns::Double(0.0),
"Linear solver iterations (multiples of the system matrix size)");
prm.declare_entry("Use static condensation", "true",
Patterns::Bool(),
"Solve the full block system or a reduced problem");
prm.declare_entry("Preconditioner type", "ssor",
Patterns::Selection("jacobi|ssor"),
"Type of preconditioner");
prm.declare_entry("Preconditioner relaxation", "0.65",
Patterns::Double(0.0),
"Preconditioner relaxation value");
}
prm.leave_subsection();
}
void LinearSolver::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Linear solver");
{
type_lin = prm.get("Solver type");
tol_lin = prm.get_double("Residual");
max_iterations_lin = prm.get_double("Max iteration multiplier");
use_static_condensation = prm.get_bool("Use static condensation");
preconditioner_type = prm.get("Preconditioner type");
preconditioner_relaxation = prm.get_double("Preconditioner relaxation");
}
prm.leave_subsection();
}
// @sect4{Nonlinear solver} ================================================================
void NonlinearSolver::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Nonlinear solver");
{
prm.declare_entry("Max iterations Newton-Raphson", "10",
Patterns::Integer(0),
"Number of Newton-Raphson iterations allowed");
prm.declare_entry("Tolerance force", "1.0e-9",
Patterns::Double(0.0),
"Force residual tolerance");
prm.declare_entry("Tolerance displacement", "1.0e-6",
Patterns::Double(0.0),
"Displacement error tolerance");
}
prm.leave_subsection();
}
void NonlinearSolver::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Nonlinear solver");
{
max_iterations_NR = prm.get_integer("Max iterations Newton-Raphson");
tol_f = prm.get_double("Tolerance force");
tol_u = prm.get_double("Tolerance displacement");
}
prm.leave_subsection();
}
// @sect4{Time} ===============================================================================
void Time::declare_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Time");
{
prm.declare_entry("End time", "1",
Patterns::Double(),
"End time");
prm.declare_entry("Time step size", "0.1",
Patterns::Double(),
"Time step size");
}
prm.leave_subsection();
}
void Time::parse_parameters(ParameterHandler &prm)
{
prm.enter_subsection("Time");
{
end_time = prm.get_double("End time");
delta_t = prm.get_double("Time step size");
}
prm.leave_subsection();
}
// @sect4{All parameters} ======================================================================
AllParameters::AllParameters(const std::string &input_file)
{
ParameterHandler prm;
declare_parameters(prm);
prm.parse_input(input_file);
parse_parameters(prm);
}
void AllParameters::declare_parameters(ParameterHandler &prm)
{
FESystem::declare_parameters(prm);
Geometry::declare_parameters(prm);
Materials::declare_parameters(prm);
LinearSolver::declare_parameters(prm);
NonlinearSolver::declare_parameters(prm);
Time::declare_parameters(prm);
}
void AllParameters::parse_parameters(ParameterHandler &prm)
{
FESystem::parse_parameters(prm);
Geometry::parse_parameters(prm);
Materials::parse_parameters(prm);
LinearSolver::parse_parameters(prm);
NonlinearSolver::parse_parameters(prm);
Time::parse_parameters(prm);
}
}