::type matrix_type;
//! \brief The domain type of the preconditioner.
typedef X domain_type;
//! \brief The range type of the preconditioner.
typedef Y range_type;
//! \brief The field type of the preconditioner.
typedef typename X::field_type field_type;
// define the category
enum {
//! \brief The category the preconditioner is part of.
category = std::is_same::value?
Dune::SolverCategory::sequential:Dune::SolverCategory::overlapping
};
//! \brief Elliptic Operator
typedef typename CPRSelector::Operator Operator;
//! \brief preconditioner for the whole system (here either ILU(0) or ILU(n)
typedef Dune::Preconditioner WholeSystemPreconditioner;
//! \brief the ilu-0 preconditioner used the for the elliptic system
typedef typename CPRSelector::EllipticPreconditioner
EllipticPreconditioner;
//! \brief type of the unique pointer to the ilu-0 preconditioner
//! used the for the elliptic system
typedef typename CPRSelector::EllipticPreconditionerPointer
EllipticPreconditionerPointer;
//! \brief amg preconditioner for the elliptic system
typedef EllipticPreconditioner Smoother;
typedef Dune::Amg::AMG AMG;
/*! \brief Constructor.
Constructor gets all parameters to operate the prec.
\param A The matrix to operate on.
\param Ae The top-left elliptic part of A.
\param relax The ILU0 relaxation factor.
\param useAMG if true, AMG is used as a preconditioner for the elliptic sub-system, otherwise ilu-0 (default)
\param useBiCG if true, BiCG solver is used (default), otherwise CG solver
\param paralleInformation The information about the parallelization, if this is a
parallel run
*/
CPRPreconditioner (const M& A, const M& Ae, const field_type relax,
const unsigned int ilu_n,
const bool useAMG,
const bool useBiCG,
const ParallelInformation& comm=ParallelInformation())
: A_(A),
Ae_(Ae),
de_( Ae_.N() ),
ve_( Ae_.M() ),
dmodified_( A_.N() ),
opAe_(CPRSelector::makeOperator(Ae_, comm)),
precond_(), // ilu0 preconditioner for elliptic system
amg_(), // amg preconditioner for elliptic system
pre_(), // copy A will be made be the preconditioner
vilu_( A_.N() ),
relax_(relax),
use_bicg_solver_( useBiCG ),
comm_(comm)
{
// create appropriate preconditioner for elliptic system
createPreconditioner( useAMG, comm );
if( ilu_n == 0 ) {
pre_ = createILU0Ptr( A_, relax_, comm );
}
else {
pre_ = createILUnPtr( A_, ilu_n, relax_, comm );
}
}
/*!
\brief Prepare the preconditioner.
\copydoc Preconditioner::pre(X&,Y&)
*/
virtual void pre (X& /*x*/, Y& /*b*/)
{
}
/*!
\brief Apply the preconditoner.
\copydoc Preconditioner::apply(X&,const Y&)
*/
virtual void apply (X& v, const Y& d)
{
// Extract part of d corresponding to elliptic part.
// Note: Assumes that the elliptic part comes first.
std::copy_n(d.begin(), de_.size(), de_.begin());
// Solve elliptic part, extend solution to full.
// reset result
ve_ = 0;
solveElliptic( ve_, de_ );
//reset return value
v = 0.0;
// Again assuming that the elliptic part comes first.
std::copy(ve_.begin(), ve_.end(), v.begin());
// Subtract elliptic residual from initial residual.
// dmodified = d - A * vfull
dmodified_ = d;
A_.mmv(v, dmodified_);
// Apply Preconditioner for whole system (relax will be applied already)
pre_->apply( vilu_, dmodified_);
// don't apply relaxation if relax_ == 1
if( std::abs( relax_ - 1.0 ) < 1e-12 ) {
v += vilu_;
}
else {
v *= relax_;
v += vilu_;
}
}
/*!
\brief Clean up.
\copydoc Preconditioner::post(X&)
*/
virtual void post (X& /*x*/)
{
}
protected:
void solveElliptic(Y& x, Y& de)
{
// Linear solver parameters
const double tolerance = 1e-4;
const int maxit = 5000;
const int verbosity = 0;
// operator result containing iterations etc.
Dune::InverseOperatorResult result;
// the scalar product chooser
typedef Dune::ScalarProductChooser
ScalarProductChooser;
// the scalar product.
std::unique_ptr
sp(ScalarProductChooser::construct(comm_));
if( amg_ )
{
// Solve system with AMG
if( use_bicg_solver_ ) {
Dune::BiCGSTABSolver linsolve(*opAe_, *sp, (*amg_), tolerance, maxit, verbosity);
linsolve.apply(x, de, result);
}
else {
Dune::CGSolver linsolve(*opAe_, *sp, (*amg_), tolerance, maxit, verbosity);
linsolve.apply(x, de, result);
}
}
else
{
assert( precond_ );
// Solve system with ILU-0
if( use_bicg_solver_ ) {
Dune::BiCGSTABSolver linsolve(*opAe_, *sp, (*precond_), tolerance, maxit, verbosity);
linsolve.apply(x, de, result);
}
else {
Dune::CGSolver linsolve(*opAe_, *sp, (*precond_), tolerance, maxit, verbosity);
linsolve.apply(x, de, result);
}
}
if (!result.converged) {
OPM_THROW(std::runtime_error, "CPRPreconditioner failed to solve elliptic subsystem.");
}
}
//! \brief The matrix for the full linear problem.
const matrix_type& A_;
//! \brief The elliptic part of the matrix.
const matrix_type& Ae_;
//! \brief temporary variables for elliptic solve
Y de_, ve_, dmodified_;
//! \brief elliptic operator
std::unique_ptr opAe_;
//! \brief ILU0 preconditioner for the elliptic system
EllipticPreconditionerPointer precond_;
//! \brief AMG preconditioner with ILU0 smoother
std::unique_ptr< AMG > amg_;
//! \brief The preconditioner for the whole system
//!
//! We have to use a shared_ptr instead of a unique_ptr
//! as we need to use a custom allocator based on dynamic
//! information. But for unique_ptr the type of this deleter
//! has to be available at coompile time.
std::shared_ptr< WholeSystemPreconditioner > pre_;
//! \brief temporary variables for ILU solve
Y vilu_;
//! \brief The relaxation factor to use.
field_type relax_;
//! \brief true if ISTL BiCGSTABSolver is used, otherwise ISTL CGSolver is used
const bool use_bicg_solver_;
//! \brief The information about the parallelization
const P& comm_;
protected:
void createPreconditioner( const bool amg, const P& comm )
{
if( amg )
{
typedef Dune::Amg::CoarsenCriterion< Dune::Amg::SymmetricCriterion > Criterion;
typedef typename Dune::Amg::SmootherTraits::Arguments SmootherArgs;
SmootherArgs smootherArgs;
smootherArgs.iterations = 1;
smootherArgs.relaxationFactor = relax_;
int coarsenTarget=1200;
Criterion criterion(15,coarsenTarget);
criterion.setDebugLevel( 0 ); // no debug information, 1 for printing hierarchy information
criterion.setDefaultValuesIsotropic(2);
criterion.setAlpha(.67);
criterion.setBeta(1.0e-6);
criterion.setMaxLevel(10);
amg_ = std::unique_ptr< AMG > (new AMG(*opAe_, criterion, smootherArgs));
}
else
{
precond_ = createEllipticPreconditionerPointer( Ae_, relax_, comm);
}
}
};
} // namespace Opm
#endif // OPM_CPRPRECONDITIONER_HEADER_INCLUDED