Merge pull request #624 from GitPaean/polymer_linearsolver

making the interleaved solver works for blackoil polymer simulator.
This commit is contained in:
Atgeirr Flø Rasmussen 2016-03-31 15:39:54 +02:00
commit 6549658622
3 changed files with 53 additions and 2 deletions

View File

@ -38,6 +38,10 @@ namespace Opm
{
protected:
using Base = FlowMainBase<FlowMainPolymer<Grid, Simulator>, Grid, Simulator>;
using Base::eclipse_state_;
using Base::param_;
using Base::fis_solver_;
using Base::parallel_information_;
friend Base;
// Set in setupGridAndProps()
@ -90,10 +94,34 @@ namespace Opm
// Setup linear solver.
// Writes to:
// fis_solver_
// Currently, the CPR solver is not ready for polymer solver yet
void setupLinearSolver()
{
OPM_MESSAGE("Caution: polymer solver currently only works with direct linear solver.");
Base::fis_solver_.reset(new NewtonIterationBlackoilSimple(Base::param_, Base::parallel_information_));
const std::string cprSolver = "cpr";
const std::string interleavedSolver = "interleaved";
const std::string directSolver = "direct";
const std::string flowDefaultSolver = interleavedSolver;
std::shared_ptr<const Opm::SimulationConfig> simCfg = eclipse_state_->getSimulationConfig();
std::string solver_approach = flowDefaultSolver;
if (param_.has("solver_approach")) {
solver_approach = param_.template get<std::string>("solver_approach");
} else {
if (simCfg->useCPR()) {
solver_approach = cprSolver;
}
}
if (solver_approach == cprSolver) {
OPM_THROW( std::runtime_error , "CPR solver is not ready for use with polymer solver yet.");
} else if (solver_approach == interleavedSolver) {
fis_solver_.reset(new NewtonIterationBlackoilInterleaved(param_, parallel_information_));
} else if (solver_approach == directSolver) {
fis_solver_.reset(new NewtonIterationBlackoilSimple(param_, parallel_information_));
} else {
OPM_THROW( std::runtime_error , "Internal error - solver approach " << solver_approach << " not recognized.");
}
}

View File

@ -234,6 +234,8 @@ namespace Opm {
const SolutionState& state,
WellState& xw);
void updateEquationsScaling();
void
computeMassFlux(const int actph ,
const V& transi,

View File

@ -107,6 +107,7 @@ namespace Opm {
if (!active_[Water]) {
OPM_THROW(std::logic_error, "Polymer must solved in water!\n");
}
residual_.matbalscale.resize(fluid_.numPhases() + 1, 1.1169); // use the same as the water phase
// If deck has polymer, residual_ should contain polymer equation.
rq_.resize(fluid_.numPhases() + 1);
residual_.material_balance_eq.resize(fluid_.numPhases() + 1, ADB::null());
@ -344,6 +345,26 @@ namespace Opm {
residual_.material_balance_eq[poly_pos_] = pvdt_ * (rq_[poly_pos_].accum[1] - rq_[poly_pos_].accum[0])
+ ops_.div*rq_[poly_pos_].mflux;
}
if (param_.update_equations_scaling_) {
updateEquationsScaling();
}
}
template <class Grid>
void BlackoilPolymerModel<Grid>::updateEquationsScaling()
{
Base::updateEquationsScaling();
if (has_polymer_) {
const int water_pos = fluid_.phaseUsage().phase_pos[Water];
residual_.matbalscale[poly_pos_] = residual_.matbalscale[water_pos];
}
}