24 #ifndef RICHARDSON_LINEAR_SOLVER_H_
25 #define RICHARDSON_LINEAR_SOLVER_H_
27 #include <deal.II/lac/vector.h>
28 #include <deal.II/lac/block_sparsity_pattern.h>
29 #include <deal.II/lac/block_sparse_matrix.h>
30 #include <deal.II/lac/compressed_simple_sparsity_pattern.h>
31 #include <deal.II/lac/solver_richardson.h>
32 #include <deal.II/lac/precondition.h>
33 #include <deal.II/lac/full_matrix.h>
35 #include <deal.II/dofs/dof_tools.h>
37 #include <deal.II/numerics/vector_tools.h>
56 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
69 template<
typename PROBLEM>
90 template<
typename PROBLEM,
typename INTEGRATOR>
91 void Solve(PROBLEM& pde,INTEGRATOR& integr, VECTOR &rhs, VECTOR &solution,
bool force_matrix_build=
false);
96 SPARSITYPATTERN sparsity_pattern_;
99 double linear_global_tol_, linear_tol_;
105 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
108 param_reader.
SetSubsection(
"richardsonwithmatrix parameters");
109 param_reader.
declare_entry(
"linear_global_tol",
"1.e-10",Patterns::Double(0),
"global tolerance for the richardson iteration");
110 param_reader.
declare_entry(
"linear_maxiter",
"1000",Patterns::Integer(0),
"maximal number of richardson steps");
114 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
118 param_reader.
SetSubsection(
"richardsonwithmatrix parameters");
119 linear_global_tol_ = param_reader.
get_double (
"linear_global_tol");
120 linear_maxiter_ = param_reader.
get_integer (
"linear_maxiter");
125 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
132 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
133 template<
typename PROBLEM>
137 pde.ComputeSparsityPattern(sparsity_pattern_);
138 matrix_.reinit(sparsity_pattern_);
143 template <
typename PRECONDITIONER,
typename SPARSITYPATTERN,
typename MATRIX,
typename VECTOR>
144 template<
typename PROBLEM,
typename INTEGRATOR>
149 bool force_matrix_build)
151 if(force_matrix_build)
153 integr.ComputeMatrix (pde,matrix_);
156 integr.ApplyNewtonBoundaryValues(pde,matrix_,rhs,solution);
158 dealii::SolverControl solver_control (linear_maxiter_, linear_global_tol_,
false,
false);
160 dealii::SolverRichardson<VECTOR> richardson(solver_control);
161 PRECONDITIONER precondition;
162 precondition.initialize(matrix_);
163 richardson.solve (matrix_, solution, rhs,
166 pde.GetDoFConstraints().distribute(solution);
void ReInit(PROBLEM &pde)
Definition: richardsonlinearsolver.h:134
double get_double(const std::string &entry_name)
Definition: parameterreader.h:115
Definition: parameterreader.h:36
void declare_entry(const std::string &entry, const std::string &default_value, const Patterns::PatternBase &pattern=Patterns::Anything(), const std::string &documentation=std::string())
Definition: parameterreader.h:98
void Solve(PROBLEM &pde, INTEGRATOR &integr, VECTOR &rhs, VECTOR &solution, bool force_matrix_build=false)
Definition: richardsonlinearsolver.h:145
RichardsonLinearSolverWithMatrix(ParameterReader ¶m_reader)
Definition: richardsonlinearsolver.h:116
~RichardsonLinearSolverWithMatrix()
Definition: richardsonlinearsolver.h:126
Definition: richardsonlinearsolver.h:57
int get_integer(const std::string &entry_name)
Definition: parameterreader.h:126
static void declare_params(ParameterReader ¶m_reader)
Definition: richardsonlinearsolver.h:106
void SetSubsection(const std::string subsection)
Definition: parameterreader.h:93