typename PreconditionerTraits::PointerType
makePreconditioner(O& opA, double relax, const C& comm, int iterations=1)
{
typedef typename SmootherChooser
::Type SmootherType;
typedef typename PreconditionerTraits
::PointerType PointerType;
typename Dune::Amg::SmootherTraits::Arguments args;
typename Dune::Amg::ConstructionTraits::Arguments cargs;
cargs.setMatrix(opA.getmat());
args.iterations=iterations;
args.relaxationFactor=relax;
cargs.setArgs(args);
cargs.setComm(comm);
return PointerType(Dune::Amg::ConstructionTraits::construct(cargs));
}
template
LinearSolverInterface::LinearSolverReport
solveCG_ILU0(O& opA, Vector& x, Vector& b, S& sp, const C& comm, double tolerance, int maxit, int verbosity)
{
// Construct preconditioner.
typedef Dune::SeqILU0 Preconditioner;
auto precond = makePreconditioner(opA, 1.0, comm);
// Construct linear solver.
Dune::CGSolver linsolve(opA, sp, *precond, tolerance, maxit, verbosity);
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
// Output results.
LinearSolverInterface::LinearSolverReport res;
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
return res;
}
#define FIRST_DIAGONAL 1
#define SYMMETRIC 1
#define SMOOTHER_ILU 0
#define ANISOTROPIC_3D 0
template
void setUpCriterion(C& criterion, double linsolver_prolongate_factor,
int verbosity, std::size_t linsolver_smooth_steps)
{
criterion.setDebugLevel(verbosity);
#if ANISOTROPIC_3D
criterion.setDefaultValuesAnisotropic(3, 2);
#endif
criterion.setProlongationDampingFactor(linsolver_prolongate_factor);
criterion.setNoPreSmoothSteps(linsolver_smooth_steps);
criterion.setNoPostSmoothSteps(linsolver_smooth_steps);
criterion.setGamma(1); // V-cycle; this is the default
}
template
LinearSolverInterface::LinearSolverReport
solveCG_AMG(O& opA, Vector& x, Vector& b, S& sp, const C& comm, double tolerance, int maxit, int verbosity,
double linsolver_prolongate_factor, int linsolver_smooth_steps)
{
// Solve with AMG solver.
#if FIRST_DIAGONAL
typedef Dune::Amg::FirstDiagonal CouplingMetric;
#else
typedef Dune::Amg::RowSum CouplingMetric;
#endif
#if SYMMETRIC
typedef Dune::Amg::SymmetricCriterion CriterionBase;
#else
typedef Dune::Amg::UnSymmetricCriterion CriterionBase;
#endif
#if SMOOTHER_ILU
typedef Dune::SeqILU0 SeqSmoother;
#else
typedef Dune::SeqSOR SeqSmoother;
#endif
typedef typename SmootherChooser::Type Smoother;
typedef Dune::Amg::CoarsenCriterion Criterion;
typedef Dune::Amg::AMG Precond;
// Construct preconditioner.
Criterion criterion;
typename Precond::SmootherArgs smootherArgs;
setUpCriterion(criterion, linsolver_prolongate_factor, verbosity,
linsolver_smooth_steps);
Precond precond(opA, criterion, smootherArgs, comm);
// Construct linear solver.
Dune::CGSolver linsolve(opA, sp, precond, tolerance, maxit, verbosity);
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
// Output results.
LinearSolverInterface::LinearSolverReport res;
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
return res;
}
#if defined(HAS_DUNE_FAST_AMG) || DUNE_VERSION_NEWER(DUNE_ISTL, 2, 3)
template
LinearSolverInterface::LinearSolverReport
solveKAMG(O& opA, Vector& x, Vector& b, S& /* sp */, const C& /* comm */, double tolerance, int maxit, int verbosity,
double linsolver_prolongate_factor, int linsolver_smooth_steps)
{
// Solve with AMG solver.
Dune::MatrixAdapter sOpA(opA.getmat());
#if FIRST_DIAGONAL
typedef Dune::Amg::FirstDiagonal CouplingMetric;
#else
typedef Dune::Amg::RowSum CouplingMetric;
#endif
#if SYMMETRIC
typedef Dune::Amg::SymmetricCriterion CriterionBase;
#else
typedef Dune::Amg::UnSymmetricCriterion CriterionBase;
#endif
#if SMOOTHER_ILU
typedef Dune::SeqILU0 Smoother;
#else
typedef Dune::SeqSOR Smoother;
#endif
typedef Dune::Amg::CoarsenCriterion Criterion;
typedef Dune::Amg::KAMG Precond;
// Construct preconditioner.
Precond::SmootherArgs smootherArgs;
Criterion criterion;
setUpCriterion(criterion, linsolver_prolongate_factor, verbosity,
linsolver_smooth_steps);
Precond precond(sOpA, criterion, smootherArgs);
// Construct linear solver.
Dune::GeneralizedPCGSolver linsolve(sOpA, precond, tolerance, maxit, verbosity);
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
// Output results.
LinearSolverInterface::LinearSolverReport res;
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
return res;
}
template
LinearSolverInterface::LinearSolverReport
solveFastAMG(O& opA, Vector& x, Vector& b, S& /* sp */, const C& /* comm */, double tolerance, int maxit, int verbosity,
double linsolver_prolongate_factor)
{
// Solve with AMG solver.
typedef Dune::MatrixAdapter Operator;
Operator sOpA(opA.getmat());
#if FIRST_DIAGONAL
typedef Dune::Amg::FirstDiagonal CouplingMetric;
#else
typedef Dune::Amg::RowSum CouplingMetric;
#endif
#if SYMMETRIC
typedef Dune::Amg::AggregationCriterion > CriterionBase;
#else
typedef Dune::Amg::AggregationCriterion > CriterionBase;
#endif
typedef Dune::Amg::CoarsenCriterion Criterion;
typedef Dune::Amg::FastAMG Precond;
// Construct preconditioner.
Criterion criterion;
const int smooth_steps = 1;
setUpCriterion(criterion, linsolver_prolongate_factor, verbosity, smooth_steps);
Dune::Amg::Parameters parms;
parms.setDebugLevel(verbosity);
parms.setNoPreSmoothSteps(smooth_steps);
parms.setNoPostSmoothSteps(smooth_steps);
parms.setProlongationDampingFactor(linsolver_prolongate_factor);
Precond precond(sOpA, criterion, parms);
// Construct linear solver.
Dune::GeneralizedPCGSolver linsolve(sOpA, precond, tolerance, maxit, verbosity);
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
// Output results.
LinearSolverInterface::LinearSolverReport res;
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
return res;
}
#endif
template
LinearSolverInterface::LinearSolverReport
solveBiCGStab_ILU0(O& opA, Vector& x, Vector& b, S& sp, const C& comm, double tolerance, int maxit, int verbosity)
{
// Construct preconditioner.
typedef Dune::SeqILU0 Preconditioner;
auto precond = makePreconditioner(opA, 1.0, comm);
// Construct linear solver.
Dune::BiCGSTABSolver linsolve(opA, sp, *precond, tolerance, maxit, verbosity);
// Solve system.
Dune::InverseOperatorResult result;
linsolve.apply(x, b, result);
// Output results.
LinearSolverInterface::LinearSolverReport res;
res.converged = result.converged;
res.iterations = result.iterations;
res.residual_reduction = result.reduction;
return res;
}
} // anonymous namespace
} // namespace Opm