Add choice of solution strategies to simulators.

The following is changed in this commit:
 - The constructor for NewtonIterationBlackoilSimple now takes
   a parameter object instead of a linear solver.
 - The fully implicit black-oil simulators can now use the CPR
   preconditioning strategy (by passing use_cpr=true) or the
   simple strategy (the default).
Note that as of this commit, the CPR preconditioning still has
not been implemented properly, and behaves just like the simple
strategy.
This commit is contained in:
Atgeirr Flø Rasmussen 2014-05-12 11:10:59 +02:00
parent 46e4ea7736
commit 54c07f3d0c
5 changed files with 31 additions and 18 deletions

View File

@ -38,6 +38,7 @@
#include <opm/core/linalg/LinearSolverFactory.hpp>
#include <opm/autodiff/NewtonIterationBlackoilSimple.hpp>
#include <opm/autodiff/NewtonIterationBlackoilCPR.hpp>
#include <opm/core/simulator/BlackoilState.hpp>
#include <opm/autodiff/WellStateFullyImplicitBlackoil.hpp>
@ -141,9 +142,13 @@ try
bool use_gravity = (gravity[0] != 0.0 || gravity[1] != 0.0 || gravity[2] != 0.0);
const double *grav = use_gravity ? &gravity[0] : 0;
// Linear solver.
LinearSolverFactory linsolver(param);
NewtonIterationBlackoilSimple fis_solver(linsolver);
// Solver for Newton iterations.
std::unique_ptr<NewtonIterationBlackoilInterface> fis_solver;
if (param.getDefault("use_cpr", false)) {
fis_solver.reset(new NewtonIterationBlackoilCPR(param));
} else {
fis_solver.reset(new NewtonIterationBlackoilSimple(param));
}
// Write parameters used for later reference.
bool output = param.getDefault("output", true);
@ -212,7 +217,7 @@ try
*new_props,
rock_comp->isActive() ? rock_comp.get() : 0,
wells,
fis_solver,
*fis_solver,
grav);
SimulatorReport episodeReport = simulator.run(simtimer, state, well_state);

View File

@ -60,6 +60,7 @@
#include <opm/core/linalg/LinearSolverFactory.hpp>
#include <opm/autodiff/NewtonIterationBlackoilSimple.hpp>
#include <opm/autodiff/NewtonIterationBlackoilCPR.hpp>
#include <opm/core/simulator/BlackoilState.hpp>
#include <opm/autodiff/WellStateFullyImplicitBlackoil.hpp>
@ -194,9 +195,13 @@ try
bool use_gravity = (gravity[0] != 0.0 || gravity[1] != 0.0 || gravity[2] != 0.0);
const double *grav = use_gravity ? &gravity[0] : 0;
// Linear solver.
LinearSolverFactory linsolver(param);
NewtonIterationBlackoilSimple fis_solver(linsolver);
// Solver for Newton iterations.
std::unique_ptr<NewtonIterationBlackoilInterface> fis_solver;
if (param.getDefault("use_cpr", false)) {
fis_solver.reset(new NewtonIterationBlackoilCPR(param));
} else {
fis_solver.reset(new NewtonIterationBlackoilSimple(param));
}
// Write parameters used for later reference.
bool output = param.getDefault("output", true);
@ -271,7 +276,7 @@ try
*new_props,
rock_comp->isActive() ? rock_comp.get() : 0,
wells,
fis_solver,
*fis_solver,
grav);
SimulatorReport episodeReport = simulator.run(simtimer, state, well_state);

View File

@ -102,8 +102,7 @@ try
double grav[] = { 0.0, 0.0 };
Opm::DerivedGeology geo(*g, props, grav);
Opm::LinearSolverFactory linsolver(param);
Opm::NewtonIterationBlackoilSimple fis_solver(linsolver);
Opm::NewtonIterationBlackoilSimple fis_solver(param);
Opm::FullyImplicitBlackoilSolver<UnstructuredGrid> solver(*g, props, geo, 0, *wells, fis_solver);

View File

@ -22,15 +22,16 @@
#include <opm/autodiff/NewtonIterationBlackoilSimple.hpp>
#include <opm/autodiff/AutoDiffHelpers.hpp>
#include <opm/core/utility/ErrorMacros.hpp>
#include <opm/core/linalg/LinearSolverFactory.hpp>
namespace Opm
{
/// Construct a system solver.
/// \param[in] linsolver linear solver to use
NewtonIterationBlackoilSimple::NewtonIterationBlackoilSimple(const LinearSolverInterface& linsolver)
: linsolver_(linsolver)
NewtonIterationBlackoilSimple::NewtonIterationBlackoilSimple(const parameter::ParameterGroup& param)
{
linsolver_.reset(new LinearSolverFactory(param));
}
/// Solve the linear system Ax = b, with A being the
@ -54,9 +55,9 @@ namespace Opm
SolutionVector dx(SolutionVector::Zero(total_residual.size()));
Opm::LinearSolverInterface::LinearSolverReport rep
= linsolver_.solve(matr.rows(), matr.nonZeros(),
matr.outerIndexPtr(), matr.innerIndexPtr(), matr.valuePtr(),
total_residual.value().data(), dx.data());
= linsolver_->solve(matr.rows(), matr.nonZeros(),
matr.outerIndexPtr(), matr.innerIndexPtr(), matr.valuePtr(),
total_residual.value().data(), dx.data());
if (!rep.converged) {
OPM_THROW(std::runtime_error,
"FullyImplicitBlackoilSolver::solveJacobianSystem(): "

View File

@ -22,7 +22,9 @@
#include <opm/autodiff/NewtonIterationBlackoilInterface.hpp>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/core/linalg/LinearSolverInterface.hpp>
#include <memory>
namespace Opm
{
@ -35,8 +37,9 @@ namespace Opm
{
public:
/// Construct a system solver.
/// \param[in] linsolver linear solver to use
NewtonIterationBlackoilSimple(const LinearSolverInterface& linsolver);
/// \param[in] param parameters controlling the behaviour and
/// choice of linear solver.
NewtonIterationBlackoilSimple(const parameter::ParameterGroup& param);
/// Solve the system of linear equations Ax = b, with A being the
/// combined derivative matrix of the residual and b
@ -46,7 +49,7 @@ namespace Opm
virtual SolutionVector computeNewtonIncrement(const LinearisedBlackoilResidual& residual) const;
private:
const LinearSolverInterface& linsolver_;
std::unique_ptr<LinearSolverInterface> linsolver_;
};
} // namespace Opm