Skip to content
Snippets Groups Projects

[WIP] Feature/solverfactory

Closed Christian Engwer requested to merge feature/solverfactory into master
3 files
+ 144
0
Compare changes
  • Side-by-side
  • Inline
Files
3
@@ -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.
Loading