Skip to content
Snippets Groups Projects
Commit 542a2f30 authored by Oliver Sander's avatar Oliver Sander Committed by sander@PCPOOL.MI.FU-BERLIN.DE
Browse files

experimental wrapper for the ISTL SeqILU0 preconditioner

[[Imported from SVN: r3133]]
parent da5987c7
No related branches found
No related tags found
No related merge requests found
#ifndef ISTL_SEQILU0_STEP_HH
#define ISTL_SEQILU0_STEP_HH
/** \file
\brief A wrapper class for the ISTL SeqILU0 implementation
*/
#include <memory>
#include <dune-solvers/iterationsteps/lineariterationstep.hh>
#include <dune/istl/preconditioners.hh>
#include <dune/istl/operators.hh>
/** \brief A wrapper class for the ISTL SeqILU0 implementation
*/
template <class MatrixType, class VectorType>
class ISTLSeqILU0Step
: public LinearIterationStep<MatrixType, VectorType>
{
typedef Dune::MatrixAdapter<MatrixType,VectorType,VectorType> Operator;
typedef Dune::SeqILU0<MatrixType,VectorType,VectorType> SeqILU0;
public:
/** \brief Constructor which initializes and sets up an algebraic hierarchy
\param smootherArgs Arguments for the smoother. See the dune-istl documentation for details
\param coarseningCriterion Arguments for the coarsening. See the dune-istl documentation for details
*/
ISTLSeqILU0Step (const MatrixType* stiffnessMatrix,
VectorType& x,
VectorType& rhs,
double relaxationFactor)
: relaxationFactor_(relaxationFactor)
{
setProblem(*stiffnessMatrix, x, rhs);
seqILU0_ = std::auto_ptr<SeqILU0>(new SeqILU0(*stiffnessMatrix, relaxationFactor));
seqILU0_->pre(*this->x_,*this->rhs_);
}
/** \brief Initialize the iteration step but don't actually build the matrix hierarchy yet */
ISTLSeqILU0Step (const MatrixType* stiffnessMatrix,
VectorType& x,
VectorType& rhs)
{
setProblem(*stiffnessMatrix, x, rhs);
}
/** \brief Sets up an algebraic hierarchy
*/
virtual void preprocess();
/** \brief Perform one iteration */
virtual void iterate();
/** \brief Get the solution */
virtual VectorType getSol();
private:
/** \brief The dune-istl sequential ILU0 preconditioner */
std::auto_ptr<SeqILU0> seqILU0_;
/** \brief Some magic relaxation factor */
double relaxationFactor_;
};
template <class MatrixType, class VectorType>
void ISTLSeqILU0Step<MatrixType,VectorType>::preprocess()
{
seqILU0_ = std::auto_ptr<SeqILU0>(new SeqILU0(*this->mat_, relaxationFactor_));
seqILU0_->pre(*this->x_,*this->rhs_);
}
template <class MatrixType, class VectorType>
void ISTLSeqILU0Step<MatrixType,VectorType>::iterate()
{
seqILU0_->apply(*this->x_,*this->rhs_);
}
template <class MatrixType, class VectorType>
VectorType ISTLSeqILU0Step<MatrixType,VectorType>::getSol()
{
return *this->x_;
}
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment