Saddle-Point
Loading...
Searching...
No Matches
step-6.h
Go to the documentation of this file.
1
15#ifndef STEP_6_H
16#define STEP_6_H
17
18#include <deal.II/base/quadrature_lib.h>
19#include <deal.II/dofs/dof_handler.h>
20#include <deal.II/fe/fe_system.h>
21#include <deal.II/grid/tria.h>
22#include <deal.II/lac/affine_constraints.h>
23#include <deal.II/lac/solver_control.h>
24#include <deal.II/lac/sparse_matrix.h>
25#include <deal.II/lac/vector.h>
26#include <deal.II/base/point.h>
27#include <deal.II/base/convergence_table.h>
28#include <deal.II/base/parameter_handler.h>
29
30using namespace dealii;
31
63template <int dim>
65{
66public:
82 {
83 public:
90 Parameters();
91
102 static void declare_parameters(ParameterHandler &prm);
103
110 void get_parameters(ParameterHandler &prm);
111
114
117 unsigned int fe_degree;
118
123
126 std::string custom_rhs_0;
127
130 std::string custom_rhs_1;
131
135 std::string custom_rhs_gamma;
136 };
137
149 SaddlePointProblem(const Parameters &parameters);
150
161 void run();
162
163private:
177 void setup_system();
178
193 void assemble_system();
194
201 void assemble_system_bak();
202
212 void solve();
213
221 void refine_grid();
222
232 void output_results(const unsigned int cycle) const;
233
243 void output_exact_results(const std::string &filename) const;
244
254 double exact_solution(const Point<dim> &p) const;
255
267 void compute_error(ConvergenceTable &mytable) const;
268
279 void write_convergence_table(ConvergenceTable &mytable) const;
280
282 const Parameters &parameters;
283
286 Triangulation<dim> triangulation;
287
290 FESystem<dim> fe;
291
293 DoFHandler<dim> dof_handler;
294
296 AffineConstraints<double> constraints;
297
299 SparseMatrix<double> system_matrix;
300
302 SparsityPattern sparsity_pattern;
303
305 Vector<double> solution;
306
308 Vector<double> system_rhs;
309};
310
311
326template <int dim>
327double coefficient(const Point<dim> &p)
328{
329 if (p.square() < 1)
330 return 1;
331 else
332 return 1;
333}
334
335
336#endif // STEP_6_H
Run-time options read from a ParameterHandler input file.
Definition step-6.h:82
unsigned int n_refinement_cycles
Number of refinement cycles (global refinements) to perform.
Definition step-6.h:113
std::string custom_rhs_gamma
Symbolic expression for the boundary forcing that appears in the boundary integral,...
Definition step-6.h:135
std::string custom_rhs_0
Symbolic expression for the right-hand side of the first component, parsed at run time by FunctionPa...
Definition step-6.h:126
void get_parameters(ParameterHandler &prm)
Read the parameter values out of a (previously parsed) ParameterHandler into this Parameters object.
Definition step-6.cc:158
unsigned int fe_degree
Polynomial degree of the scalar element used for each of the two FESystem components.
Definition step-6.h:117
static void declare_parameters(ParameterHandler &prm)
Register every Parameters field with the given ParameterHandler.
Definition step-6.cc:122
Parameters()
Default-construct the Parameters with sensible fallbacks.
Definition step-6.cc:103
double convergence_tolerance
Tolerance used by the Conjugate Gradient solver: iterations stop when the residual norm is below this...
Definition step-6.h:122
std::string custom_rhs_1
Symbolic expression for the right-hand side of the second component, parsed at run time by FunctionP...
Definition step-6.h:130
Solves the two-component saddle-point Laplace system on a disk.
Definition step-6.h:65
void run()
Top-level driver that performs the adaptive-refinement loop.
Definition step-6.cc:803
double coefficient(const Point< dim > &p)
Scalar, spatially variable coefficient .
Definition step-6.h:327