Pass elliptic matrix part to CPR preconditioner.

For now, the top left part is passed with no summing of equations or checks
for diagonal dominance.
This commit is contained in:
Atgeirr Flø Rasmussen 2014-05-15 15:49:14 +02:00
parent c10a6f0804
commit bbf11c42c6
2 changed files with 30 additions and 16 deletions

View File

@ -72,15 +72,15 @@ namespace Opm
Constructor gets all parameters to operate the prec.
\param A The matrix to operate on.
\param ne The size of the elliptic top-left part
\param Ae The top-left elliptic part of A.
\param w The ILU0 relaxation factor.
*/
CPRPreconditioner (const M& A, const int ne, const field_type w)
: ILU(A) // copy A
CPRPreconditioner (const M& A, const M& Ae, const field_type relax)
: ILU_(A), // copy A (will be overwritten by ILU decomp)
Ae_(Ae),
relax_(relax)
{
_w =w;
Dune::bilu0_decomposition(ILU);
// Ae = A(pindx, pindx);
Dune::bilu0_decomposition(ILU_);
}
/*!
@ -101,8 +101,21 @@ namespace Opm
*/
virtual void apply (X& v, const Y& d)
{
Dune::bilu_backsolve(ILU,v,d);
v *= _w;
// Extract part of d corresponding to elliptic part.
Y de(Ae_.N());
std::copy_n(d.begin(), Ae_.N(), de.begin());
// Solve elliptic part, extend solution to full.
// "ve = Ae_ \ de;"
// "vfull = ve Extend full"
// Subtract elliptic residual from initial residual.
// "dmodified = d - A*vfull"
// Apply ILU0.
Dune::bilu_backsolve(ILU_, vilu, dmodified);
// "v = vfull + vilu;"
v *= relax_;
}
/*!
@ -116,12 +129,12 @@ namespace Opm
}
private:
//! \brief The relaxation factor to use.
field_type _w;
//! \brief The ILU0 decomposition of the matrix.
matrix_type ILU;
matrix_type ILU_;
//! \brief The elliptic part of the matrix.
matrix_type Ae;
matrix_type Ae_;
//! \brief The relaxation factor to use.
field_type relax_;
};
} // namespace Opm

View File

@ -158,6 +158,10 @@ namespace Opm
// Create ISTL matrix.
Mat A = makeIstlMatrix(matr);
// Create ISTL matrix for elliptic part.
const int nc = residual.material_balance_eq[0].size();
Mat Ae = makeIstlMatrix(matr.topLeftCorner(nc, nc));
// Construct operator, scalar product and vectors needed.
typedef Dune::MatrixAdapter<Mat,Vector,Vector> Operator;
Operator opA(A);
@ -173,8 +177,7 @@ namespace Opm
// typedef Dune::SeqILU0<Mat,Vector,Vector> Preconditioner;
typedef Opm::CPRPreconditioner<Mat,Vector,Vector> Preconditioner;
const double relax = 1.0;
const int nc = residual.material_balance_eq[0].size();
Preconditioner precond(A, nc, relax);
Preconditioner precond(A, Ae, relax);
// Construct linear solver.
const double tolerance = 1e-8;
@ -182,7 +185,6 @@ namespace Opm
const int verbosity = 1;
Dune::BiCGSTABSolver<Vector> linsolve(opA, sp, precond, tolerance, maxit, verbosity);
/*
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
@ -192,7 +194,6 @@ namespace Opm
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
*/
// std::copy(x.begin(), x.end(), dx.data());