Skip to content
Snippets Groups Projects
Commit 9a6a95c6 authored by Linus Seelinger's avatar Linus Seelinger Committed by Christian Engwer
Browse files

[factory] Add preconditioner / solver constructors supporting ParameterTrees

parent cfdd639e
No related branches found
No related tags found
No related merge requests found
......@@ -12,6 +12,7 @@
#include <dune/common/simd/simd.hh>
#include <dune/common/unused.hh>
#include <dune/common/parametertree.hh>
#include "preconditioner.hh"
#include "solver.hh"
......@@ -159,6 +160,27 @@ namespace Dune {
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\brief Constructor.
\param A The matrix to operate on.
\param configuration ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning
------------------|------------
iterations | The number of iterations to perform.
relaxation | Common parameter defined [here](@ref ISTL_Factory_Common_Params).
See \ref ISTL_Factory for the ParameterTree layout and examples.
*/
SeqSSOR (const M& A, const ParameterTree& configuration)
: _A_(A)
{
_n = configuration.get<field_type>("iterations");
_w = configuration.get<field_type>("relaxation");
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\brief Prepare the preconditioner.
......@@ -249,6 +271,17 @@ namespace Dune {
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\copydoc SeqSSOR::SeqSSOR(const M&,const ParameterTree&)
*/
SeqSOR (const M& A, const ParameterTree& configuration)
: _A_(A)
{
_n = configuration.get<field_type>("iterations");
_w = configuration.get<field_type>("relaxation");
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\brief Prepare the preconditioner.
......@@ -354,6 +387,17 @@ namespace Dune {
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\copydoc SeqSSOR::SeqSSOR(const M&,const ParameterTree&)
*/
SeqGS (const M& A, const ParameterTree& configuration)
: _A_(A)
{
_n = configuration.get<field_type>("iterations");
_w = configuration.get<field_type>("relaxation");
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\brief Prepare the preconditioner.
......@@ -440,6 +484,17 @@ namespace Dune {
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\copydoc SeqSSOR::SeqSSOR(const M&,const ParameterTree&)
*/
SeqJac (const M& A, const ParameterTree& configuration)
: _A_(A)
{
_n = configuration.get<field_type>("iterations");
_w = configuration.get<field_type>("relaxation");
CheckIfDiagonalPresent<M,l>::check(_A_);
}
/*!
\brief Prepare the preconditioner.
......@@ -679,6 +734,25 @@ namespace Dune {
bilu0_decomposition(ILU);
}
/*!
\brief Constructor.
\param A The matrix to operate on.
\param configuration ParameterTree containing preconditioner parameters.
ParameterTree Key | Meaning
------------------|------------
relaxation | Common parameter defined [here](@ref ISTL_Factory_Common_Params).
See \ref ISTL_Factory for the ParameterTree layout and examples.
*/
SeqILU0 (const M& A, const ParameterTree& configuration)
: ILU(A) // copy A
{
_w = configuration.get<field_type>("relaxation");
bilu0_decomposition(ILU);
}
/*!
\brief Prepare the preconditioner.
......@@ -769,6 +843,17 @@ namespace Dune {
bilu_decomposition(A,n,ILU);
}
/*!
\copydoc SeqSSOR::SeqSSOR(const M&,const ParameterTree&)
*/
SeqILUn (const M& A, const ParameterTree& configuration)
: ILU(A.N(),A.M(),M::row_wise)
{
_n = configuration.get<field_type>("iterations");
_w = configuration.get<field_type>("relaxation");
bilu_decomposition(A,_n,ILU);
}
/*!
\brief Prepare the preconditioner.
......@@ -847,6 +932,14 @@ namespace Dune {
_w(w)
{}
/*!
\copydoc SeqILU0::SeqILU0(const M&,const ParameterTree&)
*/
Richardson (const ParameterTree& configuration)
{
_w = configuration.get<field_type>("relaxation");
}
/*!
\brief Prepare the preconditioner.
......
......@@ -11,6 +11,7 @@
#include <dune/common/shared_ptr.hh>
#include <dune/common/simd/io.hh>
#include <dune/common/simd/simd.hh>
#include <dune/common/parametertree.hh>
#include "solvertype.hh"
#include "preconditioner.hh"
......@@ -266,6 +267,36 @@ namespace Dune
DUNE_THROW(InvalidSolverCategory, "LinearOperator and ScalarProduct must have the same SolverCategory!");
}
IterativeSolver (LinearOperator<X,Y>& op, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
_op(stackobject_to_shared_ptr(op)),
_prec(stackobject_to_shared_ptr(prec)),
_sp(new SeqScalarProduct<X>),
_reduction(configuration.get<real_type>("reduction")),
_maxit(configuration.get<int>("maxit")),
_verbose(configuration.get<int>("verbose")),
_category(SolverCategory::category(*_op))
{
if(SolverCategory::category(*op) != SolverCategory::sequential)
DUNE_THROW(InvalidSolverCategory, "LinearOperator has to be sequential!");
if(SolverCategory::category(*prec) != SolverCategory::sequential)
DUNE_THROW(InvalidSolverCategory, "Preconditioner has to be sequential!");
}
IterativeSolver (LinearOperator<X,Y>& op, ScalarProduct<X>& sp, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
_op(stackobject_to_shared_ptr(op)),
_prec(prec),
_sp(sp),
_reduction(configuration.get<real_type>("reduction")),
_maxit(configuration.get<int>("maxit")),
_verbose(configuration.get<int>("verbose")),
_category(SolverCategory::category(*op))
{
if(SolverCategory::category(*op) != SolverCategory::category(*prec))
DUNE_THROW(InvalidSolverCategory, "LinearOperator and Preconditioner must have the same SolverCategory!");
if(SolverCategory::category(*op) != SolverCategory::category(*sp))
DUNE_THROW(InvalidSolverCategory, "LinearOperator and ScalarProduct must have the same SolverCategory!");
}
// #warning actually we want to have this as the default and just implement the second one
// //! \copydoc InverseOperator::apply(X&,Y&,InverseOperatorResult&)
// virtual void apply (X& x, Y& b, InverseOperatorResult& res)
......
......@@ -1060,6 +1060,16 @@ namespace Dune {
_restart(restart)
{}
RestartedGMResSolver (LinearOperator<X,Y>& op, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
IterativeSolver<X,Y>::IterativeSolver(op,prec,configuration),
_restart(configuration.get<int>("restart"))
{}
RestartedGMResSolver (LinearOperator<X,Y>& op, ScalarProduct<X>& sp, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
IterativeSolver<X,Y>::IterativeSolver(op,sp,prec,configuration),
_restart(configuration.get<int>("restart"))
{}
/*!
\brief Apply inverse operator.
......@@ -1377,6 +1387,16 @@ namespace Dune {
_restart(restart)
{}
GeneralizedPCGSolver (LinearOperator<X,Y>& op, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
IterativeSolver<X,Y>::IterativeSolver(op,prec,configuration),
_restart(configuration.get<int>("restart"))
{}
GeneralizedPCGSolver (LinearOperator<X,Y>& op, ScalarProduct<X>& sp, std::shared_ptr<Preconditioner<X,X> > prec, const ParameterTree& configuration) :
IterativeSolver<X,Y>::IterativeSolver(op,sp,prec,configuration),
_restart(configuration.get<int>("restart"))
{}
/*!
\brief Apply inverse operator.
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment