2015-05-19 14:40:38 +02:00
|
|
|
/*
|
|
|
|
|
Copyright 2013, 2015 SINTEF ICT, Applied Mathematics.
|
2015-05-21 16:29:01 +02:00
|
|
|
Copyright 2014, 2015 Dr. Blatt - HPC-Simulation-Software & Services
|
|
|
|
|
Copyright 2014, 2015 Statoil ASA.
|
2015-05-19 14:40:38 +02:00
|
|
|
Copyright 2015 NTNU
|
2015-05-21 16:29:01 +02:00
|
|
|
Copyright 2015 IRIS AS
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
|
*/
|
|
|
|
|
|
|
|
|
|
#ifndef OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
|
|
|
|
|
#define OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
|
|
|
|
|
|
|
|
|
|
#include <opm/polymer/fullyimplicit/BlackoilPolymerModel.hpp>
|
|
|
|
|
|
|
|
|
|
#include <opm/autodiff/AutoDiffBlock.hpp>
|
|
|
|
|
#include <opm/autodiff/AutoDiffHelpers.hpp>
|
|
|
|
|
#include <opm/autodiff/GridHelpers.hpp>
|
|
|
|
|
#include <opm/autodiff/BlackoilPropsAdInterface.hpp>
|
|
|
|
|
#include <opm/autodiff/GeoProps.hpp>
|
|
|
|
|
#include <opm/autodiff/WellDensitySegmented.hpp>
|
|
|
|
|
|
|
|
|
|
#include <opm/core/grid.h>
|
|
|
|
|
#include <opm/core/linalg/LinearSolverInterface.hpp>
|
|
|
|
|
#include <opm/core/linalg/ParallelIstlInformation.hpp>
|
|
|
|
|
#include <opm/core/props/rock/RockCompressibility.hpp>
|
|
|
|
|
#include <opm/core/utility/ErrorMacros.hpp>
|
|
|
|
|
#include <opm/core/utility/Exceptions.hpp>
|
|
|
|
|
#include <opm/core/utility/Units.hpp>
|
|
|
|
|
#include <opm/core/well_controls.h>
|
|
|
|
|
#include <opm/core/utility/parameters/ParameterGroup.hpp>
|
|
|
|
|
|
|
|
|
|
#include <cassert>
|
|
|
|
|
#include <cmath>
|
|
|
|
|
#include <iostream>
|
|
|
|
|
#include <iomanip>
|
|
|
|
|
#include <limits>
|
|
|
|
|
|
|
|
|
|
namespace Opm {
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2015-05-24 20:38:07 +02:00
|
|
|
namespace detail {
|
2015-05-19 14:40:38 +02:00
|
|
|
|
2015-05-24 20:38:07 +02:00
|
|
|
template <class PU>
|
|
|
|
|
int polymerPos(const PU& pu)
|
|
|
|
|
{
|
|
|
|
|
const int maxnp = Opm::BlackoilPhases::MaxNumPhases;
|
|
|
|
|
int pos = 0;
|
|
|
|
|
for (int phase = 0; phase < maxnp; ++phase) {
|
|
|
|
|
if (pu.phase_used[phase]) {
|
|
|
|
|
pos++;
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
2015-05-24 20:38:07 +02:00
|
|
|
return pos;
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
2015-05-24 20:38:07 +02:00
|
|
|
} // namespace detail
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
2015-05-24 20:38:07 +02:00
|
|
|
BlackoilPolymerModel<Grid>::BlackoilPolymerModel(const typename Base::ModelParameters& param,
|
|
|
|
|
const Grid& grid,
|
|
|
|
|
const BlackoilPropsAdInterface& fluid,
|
|
|
|
|
const DerivedGeology& geo,
|
|
|
|
|
const RockCompressibility* rock_comp_props,
|
|
|
|
|
const PolymerPropsAd& polymer_props_ad,
|
|
|
|
|
const Wells* wells,
|
|
|
|
|
const NewtonIterationBlackoilInterface& linsolver,
|
|
|
|
|
const bool has_disgas,
|
|
|
|
|
const bool has_vapoil,
|
|
|
|
|
const bool has_polymer,
|
|
|
|
|
const bool terminal_output)
|
|
|
|
|
: Base(param, grid, fluid, geo, rock_comp_props, wells, linsolver,
|
|
|
|
|
has_disgas, has_vapoil, terminal_output),
|
|
|
|
|
polymer_props_ad_(polymer_props_ad),
|
|
|
|
|
has_polymer_(has_polymer),
|
|
|
|
|
poly_pos_(detail::polymerPos(fluid.phaseUsage()))
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
if (!active_[Water]) {
|
|
|
|
|
OPM_THROW(std::logic_error, "Polymer must solved in water!\n");
|
|
|
|
|
}
|
|
|
|
|
// If deck has polymer, residual_ should contain polymer equation.
|
2015-05-24 20:38:07 +02:00
|
|
|
rq_.resize(fluid_.numPhases() + 1);
|
|
|
|
|
residual_.material_balance_eq.resize(fluid_.numPhases() + 1, ADB::null());
|
2015-05-19 14:40:38 +02:00
|
|
|
assert(poly_pos_ == fluid_.numPhases());
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void
|
|
|
|
|
BlackoilPolymerModel<Grid>::
|
2015-05-19 21:29:14 +02:00
|
|
|
prepareStep(const double dt,
|
|
|
|
|
ReservoirState& reservoir_state,
|
2015-05-24 20:38:07 +02:00
|
|
|
WellState& well_state)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-24 20:38:07 +02:00
|
|
|
Base::prepareStep(dt, reservoir_state, well_state);
|
2015-05-19 21:29:14 +02:00
|
|
|
// Initial max concentration of this time step from PolymerBlackoilState.
|
|
|
|
|
cmax_ = Eigen::Map<const V>(reservoir_state.maxconcentration().data(), Opm::AutoDiffGrid::numCells(grid_));
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
2015-05-19 21:29:14 +02:00
|
|
|
void
|
2015-05-19 14:40:38 +02:00
|
|
|
BlackoilPolymerModel<Grid>::
|
2015-05-19 21:29:14 +02:00
|
|
|
afterStep(const double /* dt */,
|
|
|
|
|
ReservoirState& reservoir_state,
|
|
|
|
|
WellState& /* well_state */)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-19 21:29:14 +02:00
|
|
|
computeCmax(reservoir_state);
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void
|
|
|
|
|
BlackoilPolymerModel<Grid>::makeConstantState(SolutionState& state) const
|
|
|
|
|
{
|
2015-05-26 01:46:34 +02:00
|
|
|
Base::makeConstantState(state);
|
2015-05-19 14:40:38 +02:00
|
|
|
state.concentration = ADB::constant(state.concentration.value());
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
2015-05-26 11:49:10 +02:00
|
|
|
std::vector<V>
|
|
|
|
|
BlackoilPolymerModel<Grid>::variableStateInitials(const ReservoirState& x,
|
|
|
|
|
const WellState& xw) const
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-26 11:49:10 +02:00
|
|
|
std::vector<V> vars0 = Base::variableStateInitials(x, xw);
|
2015-05-26 12:08:30 +02:00
|
|
|
assert(int(vars0.size()) == fluid_.numPhases() + 2);
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
// Initial polymer concentration.
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
assert (not x.concentration().empty());
|
2015-05-26 11:49:10 +02:00
|
|
|
const int nc = x.concentration().size();
|
|
|
|
|
const V c = Eigen::Map<const V>(&x.concentration()[0], nc);
|
|
|
|
|
// Concentration belongs after other reservoir vars but before well vars.
|
|
|
|
|
auto concentration_pos = vars0.begin() + fluid_.numPhases();
|
|
|
|
|
assert(concentration_pos == vars0.end() - 2);
|
|
|
|
|
vars0.insert(concentration_pos, c);
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
2015-05-26 11:49:10 +02:00
|
|
|
return vars0;
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2015-05-26 11:49:10 +02:00
|
|
|
template <class Grid>
|
|
|
|
|
std::vector<int>
|
|
|
|
|
BlackoilPolymerModel<Grid>::variableStateIndices() const
|
|
|
|
|
{
|
|
|
|
|
std::vector<int> ind = Base::variableStateIndices();
|
2015-05-26 12:08:30 +02:00
|
|
|
assert(ind.size() == 5);
|
2015-05-26 12:06:17 +02:00
|
|
|
if (has_polymer_) {
|
|
|
|
|
ind.resize(6);
|
|
|
|
|
// Concentration belongs after other reservoir vars but before well vars.
|
|
|
|
|
ind[Concentration] = fluid_.numPhases();
|
|
|
|
|
// Concentration is pushing back the well vars.
|
|
|
|
|
++ind[Qs];
|
|
|
|
|
++ind[Bhp];
|
|
|
|
|
}
|
2015-05-26 11:49:10 +02:00
|
|
|
return ind;
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2015-05-26 11:49:10 +02:00
|
|
|
template <class Grid>
|
|
|
|
|
typename BlackoilPolymerModel<Grid>::SolutionState
|
|
|
|
|
BlackoilPolymerModel<Grid>::variableStateExtractVars(const ReservoirState& x,
|
|
|
|
|
const std::vector<int>& indices,
|
|
|
|
|
std::vector<ADB>& vars) const
|
|
|
|
|
{
|
|
|
|
|
SolutionState state = Base::variableStateExtractVars(x, indices, vars);
|
2015-05-26 12:06:17 +02:00
|
|
|
if (has_polymer_) {
|
|
|
|
|
state.concentration = std::move(vars[indices[Concentration]]);
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
return state;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void
|
|
|
|
|
BlackoilPolymerModel<Grid>::computeAccum(const SolutionState& state,
|
|
|
|
|
const int aix )
|
|
|
|
|
{
|
|
|
|
|
const Opm::PhaseUsage& pu = fluid_.phaseUsage();
|
|
|
|
|
|
|
|
|
|
const ADB& press = state.pressure;
|
|
|
|
|
const ADB& temp = state.temperature;
|
|
|
|
|
const std::vector<ADB>& sat = state.saturation;
|
|
|
|
|
const ADB& rs = state.rs;
|
|
|
|
|
const ADB& rv = state.rv;
|
|
|
|
|
const ADB& c = state.concentration;
|
|
|
|
|
|
|
|
|
|
const std::vector<PhasePresence> cond = phaseCondition();
|
|
|
|
|
|
|
|
|
|
const ADB pv_mult = poroMult(press);
|
|
|
|
|
|
|
|
|
|
const int maxnp = Opm::BlackoilPhases::MaxNumPhases;
|
|
|
|
|
for (int phase = 0; phase < maxnp; ++phase) {
|
|
|
|
|
if (active_[ phase ]) {
|
|
|
|
|
const int pos = pu.phase_pos[ phase ];
|
|
|
|
|
rq_[pos].b = fluidReciprocFVF(phase, state.canonical_phase_pressures[phase], temp, rs, rv, cond, cells_);
|
|
|
|
|
rq_[pos].accum[aix] = pv_mult * rq_[pos].b * sat[pos];
|
2015-05-21 16:29:01 +02:00
|
|
|
// OPM_AD_DUMP(rq_[pos].b);
|
|
|
|
|
// OPM_AD_DUMP(rq_[pos].accum[aix]);
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if (active_[ Oil ] && active_[ Gas ]) {
|
|
|
|
|
// Account for gas dissolved in oil and vaporized oil
|
|
|
|
|
const int po = pu.phase_pos[ Oil ];
|
|
|
|
|
const int pg = pu.phase_pos[ Gas ];
|
|
|
|
|
|
|
|
|
|
// Temporary copy to avoid contribution of dissolved gas in the vaporized oil
|
|
|
|
|
// when both dissolved gas and vaporized oil are present.
|
|
|
|
|
const ADB accum_gas_copy =rq_[pg].accum[aix];
|
|
|
|
|
|
|
|
|
|
rq_[pg].accum[aix] += state.rs * rq_[po].accum[aix];
|
|
|
|
|
rq_[po].accum[aix] += state.rv * accum_gas_copy;
|
2015-05-21 16:29:01 +02:00
|
|
|
// OPM_AD_DUMP(rq_[pg].accum[aix]);
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
// compute polymer properties.
|
|
|
|
|
const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
|
|
|
|
|
const ADB ads = polymer_props_ad_.adsorption(state.concentration, cmax);
|
|
|
|
|
const double rho_rock = polymer_props_ad_.rockDensity();
|
|
|
|
|
const V phi = Eigen::Map<const V>(& fluid_.porosity()[0], AutoDiffGrid::numCells(grid_), 1);
|
|
|
|
|
const double dead_pore_vol = polymer_props_ad_.deadPoreVol();
|
|
|
|
|
// compute total phases and determin polymer position.
|
|
|
|
|
rq_[poly_pos_].accum[aix] = pv_mult * rq_[pu.phase_pos[Water]].b * sat[pu.phase_pos[Water]] * c * (1. - dead_pore_vol)
|
|
|
|
|
+ pv_mult * rho_rock * (1. - phi) / phi * ads;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
2015-05-24 20:38:07 +02:00
|
|
|
void BlackoilPolymerModel<Grid>::computeCmax(ReservoirState& state)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
|
|
|
|
const int nc = AutoDiffGrid::numCells(grid_);
|
|
|
|
|
V tmp = V::Zero(nc);
|
|
|
|
|
for (int i = 0; i < nc; ++i) {
|
|
|
|
|
tmp[i] = std::max(state.maxconcentration()[i], state.concentration()[i]);
|
|
|
|
|
}
|
|
|
|
|
std::copy(&tmp[0], &tmp[0] + nc, state.maxconcentration().begin());
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void
|
|
|
|
|
BlackoilPolymerModel<Grid>::
|
2015-05-26 01:46:34 +02:00
|
|
|
assembleMassBalanceEq(const SolutionState& state)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-26 01:46:34 +02:00
|
|
|
Base::assembleMassBalanceEq(state);
|
2015-05-19 14:40:38 +02:00
|
|
|
// Add polymer equation.
|
|
|
|
|
if (has_polymer_) {
|
2015-05-19 21:29:14 +02:00
|
|
|
residual_.material_balance_eq[poly_pos_] = pvdt_ * (rq_[poly_pos_].accum[1] - rq_[poly_pos_].accum[0])
|
2015-05-19 14:40:38 +02:00
|
|
|
+ ops_.div*rq_[poly_pos_].mflux;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void BlackoilPolymerModel<Grid>::addWellEq(const SolutionState& state,
|
2015-05-24 20:38:07 +02:00
|
|
|
WellState& xw,
|
2015-05-19 21:29:14 +02:00
|
|
|
V& aliveWells)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
|
|
|
|
if( ! wellsActive() ) return ;
|
|
|
|
|
|
|
|
|
|
const int nc = Opm::AutoDiffGrid::numCells(grid_);
|
|
|
|
|
const int np = wells().number_of_phases;
|
|
|
|
|
const int nw = wells().number_of_wells;
|
|
|
|
|
const int nperf = wells().well_connpos[nw];
|
|
|
|
|
const Opm::PhaseUsage& pu = fluid_.phaseUsage();
|
|
|
|
|
V Tw = Eigen::Map<const V>(wells().WI, nperf);
|
|
|
|
|
const std::vector<int> well_cells(wells().well_cells, wells().well_cells + nperf);
|
|
|
|
|
|
|
|
|
|
// pressure diffs computed already (once per step, not changing per iteration)
|
|
|
|
|
const V& cdp = well_perforation_pressure_diffs_;
|
|
|
|
|
|
|
|
|
|
// Extract needed quantities for the perforation cells
|
|
|
|
|
const ADB& p_perfcells = subset(state.pressure, well_cells);
|
|
|
|
|
const ADB& rv_perfcells = subset(state.rv,well_cells);
|
|
|
|
|
const ADB& rs_perfcells = subset(state.rs,well_cells);
|
|
|
|
|
std::vector<ADB> mob_perfcells(np, ADB::null());
|
|
|
|
|
std::vector<ADB> b_perfcells(np, ADB::null());
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
mob_perfcells[phase] = subset(rq_[phase].mob,well_cells);
|
|
|
|
|
b_perfcells[phase] = subset(rq_[phase].b,well_cells);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Perforation pressure
|
|
|
|
|
const ADB perfpressure = (wops_.w2p * state.bhp) + cdp;
|
|
|
|
|
std::vector<double> perfpressure_d(perfpressure.value().data(), perfpressure.value().data() + nperf);
|
|
|
|
|
xw.perfPress() = perfpressure_d;
|
|
|
|
|
|
|
|
|
|
// Pressure drawdown (also used to determine direction of flow)
|
|
|
|
|
const ADB drawdown = p_perfcells - perfpressure;
|
|
|
|
|
|
|
|
|
|
// Compute vectors with zero and ones that
|
|
|
|
|
// selects the wanted quantities.
|
|
|
|
|
|
|
|
|
|
// selects injection perforations
|
|
|
|
|
V selectInjectingPerforations = V::Zero(nperf);
|
|
|
|
|
// selects producing perforations
|
|
|
|
|
V selectProducingPerforations = V::Zero(nperf);
|
|
|
|
|
for (int c = 0; c < nperf; ++c){
|
|
|
|
|
if (drawdown.value()[c] < 0)
|
|
|
|
|
selectInjectingPerforations[c] = 1;
|
|
|
|
|
else
|
|
|
|
|
selectProducingPerforations[c] = 1;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// HANDLE FLOW INTO WELLBORE
|
|
|
|
|
|
|
|
|
|
// compute phase volumetric rates at standard conditions
|
|
|
|
|
std::vector<ADB> cq_ps(np, ADB::null());
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
const ADB cq_p = -(selectProducingPerforations * Tw) * (mob_perfcells[phase] * drawdown);
|
|
|
|
|
cq_ps[phase] = b_perfcells[phase] * cq_p;
|
|
|
|
|
}
|
|
|
|
|
if (active_[Oil] && active_[Gas]) {
|
|
|
|
|
const int oilpos = pu.phase_pos[Oil];
|
|
|
|
|
const int gaspos = pu.phase_pos[Gas];
|
|
|
|
|
const ADB cq_psOil = cq_ps[oilpos];
|
|
|
|
|
const ADB cq_psGas = cq_ps[gaspos];
|
|
|
|
|
cq_ps[gaspos] += rs_perfcells * cq_psOil;
|
|
|
|
|
cq_ps[oilpos] += rv_perfcells * cq_psGas;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// HANDLE FLOW OUT FROM WELLBORE
|
|
|
|
|
|
|
|
|
|
// Using total mobilities
|
|
|
|
|
ADB total_mob = mob_perfcells[0];
|
|
|
|
|
for (int phase = 1; phase < np; ++phase) {
|
|
|
|
|
total_mob += mob_perfcells[phase];
|
|
|
|
|
}
|
|
|
|
|
// injection perforations total volume rates
|
|
|
|
|
const ADB cqt_i = -(selectInjectingPerforations * Tw) * (total_mob * drawdown);
|
|
|
|
|
|
|
|
|
|
// compute wellbore mixture for injecting perforations
|
|
|
|
|
// The wellbore mixture depends on the inflow from the reservoar
|
|
|
|
|
// and the well injection rates.
|
|
|
|
|
|
|
|
|
|
// compute avg. and total wellbore phase volumetric rates at standard conds
|
|
|
|
|
const DataBlock compi = Eigen::Map<const DataBlock>(wells().comp_frac, nw, np);
|
|
|
|
|
std::vector<ADB> wbq(np, ADB::null());
|
|
|
|
|
ADB wbqt = ADB::constant(V::Zero(nw));
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
const ADB& q_ps = wops_.p2w * cq_ps[phase];
|
|
|
|
|
const ADB& q_s = subset(state.qs, Span(nw, 1, phase*nw));
|
|
|
|
|
Selector<double> injectingPhase_selector(q_s.value(), Selector<double>::GreaterZero);
|
|
|
|
|
const int pos = pu.phase_pos[phase];
|
|
|
|
|
wbq[phase] = (compi.col(pos) * injectingPhase_selector.select(q_s,ADB::constant(V::Zero(nw)))) - q_ps;
|
|
|
|
|
wbqt += wbq[phase];
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// compute wellbore mixture at standard conditions.
|
|
|
|
|
Selector<double> notDeadWells_selector(wbqt.value(), Selector<double>::Zero);
|
|
|
|
|
std::vector<ADB> cmix_s(np, ADB::null());
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
const int pos = pu.phase_pos[phase];
|
|
|
|
|
cmix_s[phase] = wops_.w2p * notDeadWells_selector.select(ADB::constant(compi.col(pos)), wbq[phase]/wbqt);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// compute volume ratio between connection at standard conditions
|
|
|
|
|
ADB volumeRatio = ADB::constant(V::Zero(nperf));
|
|
|
|
|
const ADB d = V::Constant(nperf,1.0) - rv_perfcells * rs_perfcells;
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
ADB tmp = cmix_s[phase];
|
|
|
|
|
|
|
|
|
|
if (phase == Oil && active_[Gas]) {
|
|
|
|
|
const int gaspos = pu.phase_pos[Gas];
|
|
|
|
|
tmp = tmp - rv_perfcells * cmix_s[gaspos] / d;
|
|
|
|
|
}
|
|
|
|
|
if (phase == Gas && active_[Oil]) {
|
|
|
|
|
const int oilpos = pu.phase_pos[Oil];
|
|
|
|
|
tmp = tmp - rs_perfcells * cmix_s[oilpos] / d;
|
|
|
|
|
}
|
|
|
|
|
volumeRatio += tmp / b_perfcells[phase];
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// injecting connections total volumerates at standard conditions
|
|
|
|
|
ADB cqt_is = cqt_i/volumeRatio;
|
|
|
|
|
|
|
|
|
|
// connection phase volumerates at standard conditions
|
|
|
|
|
std::vector<ADB> cq_s(np, ADB::null());
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
cq_s[phase] = cq_ps[phase] + cmix_s[phase]*cqt_is;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Add well contributions to mass balance equations
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
residual_.material_balance_eq[phase] -= superset(cq_s[phase],well_cells,nc);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Add well contributions to polymer mass balance equation
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
const ADB mc = computeMc(state);
|
2015-05-19 21:29:14 +02:00
|
|
|
const V polyin = Eigen::Map<const V>(xw.polymerInflow().data(), nc);
|
2015-05-19 14:40:38 +02:00
|
|
|
const V poly_in_perf = subset(polyin, well_cells);
|
|
|
|
|
const V poly_mc_perf = subset(mc, well_cells).value();
|
|
|
|
|
residual_.material_balance_eq[poly_pos_] -= superset(cq_ps[pu.phase_pos[Water]] * poly_mc_perf
|
|
|
|
|
+ cmix_s[pu.phase_pos[Water]] * cqt_is*poly_in_perf,
|
|
|
|
|
well_cells, nc);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// WELL EQUATIONS
|
|
|
|
|
ADB qs = state.qs;
|
|
|
|
|
for (int phase = 0; phase < np; ++phase) {
|
|
|
|
|
qs -= superset(wops_.p2w * cq_s[phase], Span(nw, 1, phase*nw), nw*np);
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// check for dead wells (used in the well controll equations)
|
|
|
|
|
aliveWells = V::Constant(nw, 1.0);
|
|
|
|
|
for (int w = 0; w < nw; ++w) {
|
|
|
|
|
if (wbqt.value()[w] == 0) {
|
|
|
|
|
aliveWells[w] = 0.0;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Update the perforation phase rates (used to calculate the pressure drop in the wellbore)
|
|
|
|
|
V cq = superset(cq_s[0].value(), Span(nperf, np, 0), nperf*np);
|
|
|
|
|
for (int phase = 1; phase < np; ++phase) {
|
|
|
|
|
cq += superset(cq_s[phase].value(), Span(nperf, np, phase), nperf*np);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
std::vector<double> cq_d(cq.data(), cq.data() + nperf*np);
|
|
|
|
|
xw.perfPhaseRates() = cq_d;
|
|
|
|
|
|
|
|
|
|
residual_.well_flux_eq = qs;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void BlackoilPolymerModel<Grid>::updateState(const V& dx,
|
2015-05-24 20:38:07 +02:00
|
|
|
ReservoirState& reservoir_state,
|
|
|
|
|
WellState& well_state)
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
|
|
|
|
using namespace Opm::AutoDiffGrid;
|
|
|
|
|
const int np = fluid_.numPhases();
|
|
|
|
|
const int nc = numCells(grid_);
|
|
|
|
|
const int nw = wellsActive() ? wells().number_of_wells : 0;
|
|
|
|
|
const V null;
|
|
|
|
|
assert(null.size() == 0);
|
|
|
|
|
const V zero = V::Zero(nc);
|
|
|
|
|
|
|
|
|
|
// store cell status in vectors
|
|
|
|
|
V isRs = V::Zero(nc,1);
|
|
|
|
|
V isRv = V::Zero(nc,1);
|
|
|
|
|
V isSg = V::Zero(nc,1);
|
|
|
|
|
if (active_[Gas]) {
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
switch (primalVariable_[c]) {
|
|
|
|
|
case PrimalVariables::RS:
|
|
|
|
|
isRs[c] = 1;
|
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
case PrimalVariables::RV:
|
|
|
|
|
isRv[c] = 1;
|
|
|
|
|
break;
|
|
|
|
|
|
|
|
|
|
default:
|
|
|
|
|
isSg[c] = 1;
|
|
|
|
|
break;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Extract parts of dx corresponding to each part.
|
|
|
|
|
const V dp = subset(dx, Span(nc));
|
|
|
|
|
int varstart = nc;
|
|
|
|
|
const V dsw = active_[Water] ? subset(dx, Span(nc, 1, varstart)) : null;
|
|
|
|
|
varstart += dsw.size();
|
|
|
|
|
|
|
|
|
|
const V dxvar = active_[Gas] ? subset(dx, Span(nc, 1, varstart)): null;
|
|
|
|
|
varstart += dxvar.size();
|
|
|
|
|
|
|
|
|
|
const V dc = (has_polymer_) ? subset(dx, Span(nc, 1, varstart)) : null;
|
|
|
|
|
varstart += dc.size();
|
|
|
|
|
|
|
|
|
|
const V dqs = subset(dx, Span(np*nw, 1, varstart));
|
|
|
|
|
varstart += dqs.size();
|
|
|
|
|
const V dbhp = subset(dx, Span(nw, 1, varstart));
|
|
|
|
|
varstart += dbhp.size();
|
|
|
|
|
assert(varstart == dx.size());
|
|
|
|
|
|
|
|
|
|
// Pressure update.
|
|
|
|
|
const double dpmaxrel = dpMaxRel();
|
2015-05-20 08:48:54 +02:00
|
|
|
const V p_old = Eigen::Map<const V>(&reservoir_state.pressure()[0], nc, 1);
|
2015-05-19 14:40:38 +02:00
|
|
|
const V absdpmax = dpmaxrel*p_old.abs();
|
|
|
|
|
const V dp_limited = sign(dp) * dp.abs().min(absdpmax);
|
|
|
|
|
const V p = (p_old - dp_limited).max(zero);
|
2015-05-20 08:48:54 +02:00
|
|
|
std::copy(&p[0], &p[0] + nc, reservoir_state.pressure().begin());
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
// Saturation updates.
|
|
|
|
|
const Opm::PhaseUsage& pu = fluid_.phaseUsage();
|
2015-05-20 08:48:54 +02:00
|
|
|
const DataBlock s_old = Eigen::Map<const DataBlock>(& reservoir_state.saturation()[0], nc, np);
|
2015-05-19 14:40:38 +02:00
|
|
|
const double dsmax = dsMax();
|
|
|
|
|
|
|
|
|
|
V so;
|
|
|
|
|
V sw;
|
|
|
|
|
V sg;
|
|
|
|
|
|
|
|
|
|
{
|
|
|
|
|
V maxVal = zero;
|
|
|
|
|
V dso = zero;
|
|
|
|
|
if (active_[Water]){
|
|
|
|
|
maxVal = dsw.abs().max(maxVal);
|
|
|
|
|
dso = dso - dsw;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
V dsg;
|
|
|
|
|
if (active_[Gas]){
|
|
|
|
|
dsg = isSg * dxvar - isRv * dsw;
|
|
|
|
|
maxVal = dsg.abs().max(maxVal);
|
|
|
|
|
dso = dso - dsg;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
maxVal = dso.abs().max(maxVal);
|
|
|
|
|
|
|
|
|
|
V step = dsmax/maxVal;
|
|
|
|
|
step = step.min(1.);
|
|
|
|
|
|
|
|
|
|
if (active_[Water]) {
|
|
|
|
|
const int pos = pu.phase_pos[ Water ];
|
|
|
|
|
const V sw_old = s_old.col(pos);
|
|
|
|
|
sw = sw_old - step * dsw;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if (active_[Gas]) {
|
|
|
|
|
const int pos = pu.phase_pos[ Gas ];
|
|
|
|
|
const V sg_old = s_old.col(pos);
|
|
|
|
|
sg = sg_old - step * dsg;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
const int pos = pu.phase_pos[ Oil ];
|
|
|
|
|
const V so_old = s_old.col(pos);
|
|
|
|
|
so = so_old - step * dso;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Appleyard chop process.
|
|
|
|
|
auto ixg = sg < 0;
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
if (ixg[c]) {
|
|
|
|
|
sw[c] = sw[c] / (1-sg[c]);
|
|
|
|
|
so[c] = so[c] / (1-sg[c]);
|
|
|
|
|
sg[c] = 0;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
auto ixo = so < 0;
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
if (ixo[c]) {
|
|
|
|
|
sw[c] = sw[c] / (1-so[c]);
|
|
|
|
|
sg[c] = sg[c] / (1-so[c]);
|
|
|
|
|
so[c] = 0;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
auto ixw = sw < 0;
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
if (ixw[c]) {
|
|
|
|
|
so[c] = so[c] / (1-sw[c]);
|
|
|
|
|
sg[c] = sg[c] / (1-sw[c]);
|
|
|
|
|
sw[c] = 0;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
const V sumSat = sw + so + sg;
|
|
|
|
|
sw = sw / sumSat;
|
|
|
|
|
so = so / sumSat;
|
|
|
|
|
sg = sg / sumSat;
|
|
|
|
|
|
2015-05-20 08:48:54 +02:00
|
|
|
// Update the reservoir_state
|
2015-05-19 14:40:38 +02:00
|
|
|
for (int c = 0; c < nc; ++c) {
|
2015-05-20 08:48:54 +02:00
|
|
|
reservoir_state.saturation()[c*np + pu.phase_pos[ Water ]] = sw[c];
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
2015-05-20 08:48:54 +02:00
|
|
|
reservoir_state.saturation()[c*np + pu.phase_pos[ Gas ]] = sg[c];
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if (active_[ Oil ]) {
|
|
|
|
|
const int pos = pu.phase_pos[ Oil ];
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
2015-05-20 08:48:54 +02:00
|
|
|
reservoir_state.saturation()[c*np + pos] = so[c];
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Update rs and rv
|
|
|
|
|
const double drmaxrel = drMaxRel();
|
|
|
|
|
V rs;
|
|
|
|
|
if (has_disgas_) {
|
2015-05-20 08:48:54 +02:00
|
|
|
const V rs_old = Eigen::Map<const V>(&reservoir_state.gasoilratio()[0], nc);
|
2015-05-19 14:40:38 +02:00
|
|
|
const V drs = isRs * dxvar;
|
|
|
|
|
const V drs_limited = sign(drs) * drs.abs().min(rs_old.abs()*drmaxrel);
|
|
|
|
|
rs = rs_old - drs_limited;
|
|
|
|
|
}
|
|
|
|
|
V rv;
|
|
|
|
|
if (has_vapoil_) {
|
2015-05-20 08:48:54 +02:00
|
|
|
const V rv_old = Eigen::Map<const V>(&reservoir_state.rv()[0], nc);
|
2015-05-19 14:40:38 +02:00
|
|
|
const V drv = isRv * dxvar;
|
|
|
|
|
const V drv_limited = sign(drv) * drv.abs().min(rv_old.abs()*drmaxrel);
|
|
|
|
|
rv = rv_old - drv_limited;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
// Sg is used as primal variable for water only cells.
|
|
|
|
|
const double epsilon = std::sqrt(std::numeric_limits<double>::epsilon());
|
|
|
|
|
auto watOnly = sw > (1 - epsilon);
|
|
|
|
|
|
|
|
|
|
// phase translation sg <-> rs
|
|
|
|
|
std::fill(primalVariable_.begin(), primalVariable_.end(), PrimalVariables::Sg);
|
|
|
|
|
|
|
|
|
|
if (has_disgas_) {
|
|
|
|
|
const V rsSat0 = fluidRsSat(p_old, s_old.col(pu.phase_pos[Oil]), cells_);
|
|
|
|
|
const V rsSat = fluidRsSat(p, so, cells_);
|
|
|
|
|
// The obvious case
|
|
|
|
|
auto hasGas = (sg > 0 && isRs == 0);
|
|
|
|
|
|
|
|
|
|
// Set oil saturated if previous rs is sufficiently large
|
2015-05-20 08:48:54 +02:00
|
|
|
const V rs_old = Eigen::Map<const V>(&reservoir_state.gasoilratio()[0], nc);
|
2015-05-19 14:40:38 +02:00
|
|
|
auto gasVaporized = ( (rs > rsSat * (1+epsilon) && isRs == 1 ) && (rs_old > rsSat0 * (1-epsilon)) );
|
|
|
|
|
auto useSg = watOnly || hasGas || gasVaporized;
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
if (useSg[c]) {
|
|
|
|
|
rs[c] = rsSat[c];
|
|
|
|
|
} else {
|
|
|
|
|
primalVariable_[c] = PrimalVariables::RS;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// phase transitions so <-> rv
|
|
|
|
|
if (has_vapoil_) {
|
|
|
|
|
|
|
|
|
|
// The gas pressure is needed for the rvSat calculations
|
|
|
|
|
const V gaspress_old = computeGasPressure(p_old, s_old.col(Water), s_old.col(Oil), s_old.col(Gas));
|
|
|
|
|
const V gaspress = computeGasPressure(p, sw, so, sg);
|
|
|
|
|
const V rvSat0 = fluidRvSat(gaspress_old, s_old.col(pu.phase_pos[Oil]), cells_);
|
|
|
|
|
const V rvSat = fluidRvSat(gaspress, so, cells_);
|
|
|
|
|
|
|
|
|
|
// The obvious case
|
|
|
|
|
auto hasOil = (so > 0 && isRv == 0);
|
|
|
|
|
|
|
|
|
|
// Set oil saturated if previous rv is sufficiently large
|
2015-05-20 08:48:54 +02:00
|
|
|
const V rv_old = Eigen::Map<const V>(&reservoir_state.rv()[0], nc);
|
2015-05-19 14:40:38 +02:00
|
|
|
auto oilCondensed = ( (rv > rvSat * (1+epsilon) && isRv == 1) && (rv_old > rvSat0 * (1-epsilon)) );
|
|
|
|
|
auto useSg = watOnly || hasOil || oilCondensed;
|
|
|
|
|
for (int c = 0; c < nc; ++c) {
|
|
|
|
|
if (useSg[c]) {
|
|
|
|
|
rv[c] = rvSat[c];
|
|
|
|
|
} else {
|
|
|
|
|
primalVariable_[c] = PrimalVariables::RV;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
|
|
2015-05-20 08:48:54 +02:00
|
|
|
// Update the reservoir_state
|
2015-05-19 14:40:38 +02:00
|
|
|
if (has_disgas_) {
|
2015-05-20 08:48:54 +02:00
|
|
|
std::copy(&rs[0], &rs[0] + nc, reservoir_state.gasoilratio().begin());
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if (has_vapoil_) {
|
2015-05-20 08:48:54 +02:00
|
|
|
std::copy(&rv[0], &rv[0] + nc, reservoir_state.rv().begin());
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//Polymer concentration updates.
|
|
|
|
|
if (has_polymer_) {
|
2015-05-20 08:48:54 +02:00
|
|
|
const V c_old = Eigen::Map<const V>(&reservoir_state.concentration()[0], nc, 1);
|
2015-05-19 14:40:38 +02:00
|
|
|
const V c = (c_old - dc).max(zero);
|
2015-05-20 08:48:54 +02:00
|
|
|
std::copy(&c[0], &c[0] + nc, reservoir_state.concentration().begin());
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if( wellsActive() )
|
|
|
|
|
{
|
|
|
|
|
// Qs update.
|
|
|
|
|
// Since we need to update the wellrates, that are ordered by wells,
|
|
|
|
|
// from dqs which are ordered by phase, the simplest is to compute
|
|
|
|
|
// dwr, which is the data from dqs but ordered by wells.
|
|
|
|
|
const DataBlock wwr = Eigen::Map<const DataBlock>(dqs.data(), np, nw).transpose();
|
|
|
|
|
const V dwr = Eigen::Map<const V>(wwr.data(), nw*np);
|
|
|
|
|
const V wr_old = Eigen::Map<const V>(&well_state.wellRates()[0], nw*np);
|
|
|
|
|
const V wr = wr_old - dwr;
|
|
|
|
|
std::copy(&wr[0], &wr[0] + wr.size(), well_state.wellRates().begin());
|
|
|
|
|
|
|
|
|
|
// Bhp update.
|
|
|
|
|
const V bhp_old = Eigen::Map<const V>(&well_state.bhp()[0], nw, 1);
|
|
|
|
|
const V dbhp_limited = sign(dbhp) * dbhp.abs().min(bhp_old.abs()*dpmaxrel);
|
|
|
|
|
const V bhp = bhp_old - dbhp_limited;
|
|
|
|
|
std::copy(&bhp[0], &bhp[0] + bhp.size(), well_state.bhp().begin());
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// Update phase conditions used for property calculations.
|
|
|
|
|
updatePhaseCondFromPrimalVariable();
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
void
|
|
|
|
|
BlackoilPolymerModel<Grid>::computeMassFlux(const int actph ,
|
|
|
|
|
const V& transi,
|
|
|
|
|
const ADB& kr ,
|
|
|
|
|
const ADB& phasePressure,
|
|
|
|
|
const SolutionState& state)
|
|
|
|
|
{
|
|
|
|
|
const int canonicalPhaseIdx = canph_[ actph ];
|
|
|
|
|
|
|
|
|
|
const std::vector<PhasePresence> cond = phaseCondition();
|
|
|
|
|
|
|
|
|
|
const ADB tr_mult = transMult(state.pressure);
|
|
|
|
|
const ADB mu = fluidViscosity(canonicalPhaseIdx, phasePressure, state.temperature, state.rs, state.rv,cond, cells_);
|
|
|
|
|
|
|
|
|
|
rq_[ actph ].mob = tr_mult * kr / mu;
|
|
|
|
|
|
|
|
|
|
const ADB rho = fluidDensity(canonicalPhaseIdx, phasePressure, state.temperature, state.rs, state.rv,cond, cells_);
|
|
|
|
|
|
|
|
|
|
ADB& head = rq_[ actph ].head;
|
|
|
|
|
|
|
|
|
|
// compute gravity potensial using the face average as in eclipse and MRST
|
|
|
|
|
const ADB rhoavg = ops_.caver * rho;
|
|
|
|
|
|
|
|
|
|
ADB dp = ops_.ngrad * phasePressure - geo_.gravity()[2] * (rhoavg * (ops_.ngrad * geo_.z().matrix()));
|
|
|
|
|
|
|
|
|
|
if (use_threshold_pressure_) {
|
|
|
|
|
applyThresholdPressures(dp);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
head = transi*dp;
|
|
|
|
|
|
|
|
|
|
if (canonicalPhaseIdx == Water) {
|
|
|
|
|
if(has_polymer_) {
|
|
|
|
|
const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
|
|
|
|
|
const ADB mc = computeMc(state);
|
|
|
|
|
ADB krw_eff = polymer_props_ad_.effectiveRelPerm(state.concentration,
|
|
|
|
|
cmax,
|
|
|
|
|
kr,
|
|
|
|
|
state.saturation[canonicalPhaseIdx]);
|
|
|
|
|
ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(state.concentration, mu.value().data());
|
|
|
|
|
rq_[actph].mob = tr_mult * krw_eff * inv_wat_eff_visc;
|
|
|
|
|
rq_[poly_pos_].mob = tr_mult * mc * krw_eff * inv_wat_eff_visc;
|
|
|
|
|
rq_[poly_pos_].b = rq_[actph].b;
|
|
|
|
|
rq_[poly_pos_].head = rq_[actph].head;
|
|
|
|
|
UpwindSelector<double> upwind(grid_, ops_, rq_[poly_pos_].head.value());
|
|
|
|
|
rq_[poly_pos_].mflux = upwind.select(rq_[poly_pos_].b * rq_[poly_pos_].mob) * rq_[poly_pos_].head;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
//head = transi*(ops_.ngrad * phasePressure) + gflux;
|
|
|
|
|
|
|
|
|
|
UpwindSelector<double> upwind(grid_, ops_, head.value());
|
|
|
|
|
|
|
|
|
|
const ADB& b = rq_[ actph ].b;
|
|
|
|
|
const ADB& mob = rq_[ actph ].mob;
|
|
|
|
|
rq_[ actph ].mflux = upwind.select(b * mob) * head;
|
2015-05-21 16:29:01 +02:00
|
|
|
// OPM_AD_DUMP(rq_[ actph ].mob);
|
|
|
|
|
// OPM_AD_DUMP(rq_[ actph ].mflux);
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
std::vector<double>
|
|
|
|
|
BlackoilPolymerModel<Grid>::computeResidualNorms() const
|
|
|
|
|
{
|
|
|
|
|
std::vector<double> residualNorms;
|
|
|
|
|
|
|
|
|
|
std::vector<ADB>::const_iterator massBalanceIt = residual_.material_balance_eq.begin();
|
|
|
|
|
const std::vector<ADB>::const_iterator endMassBalanceIt = residual_.material_balance_eq.end();
|
|
|
|
|
|
|
|
|
|
for (; massBalanceIt != endMassBalanceIt; ++massBalanceIt) {
|
|
|
|
|
const double massBalanceResid = detail::infinityNorm( (*massBalanceIt) );
|
|
|
|
|
if (!std::isfinite(massBalanceResid)) {
|
|
|
|
|
OPM_THROW(Opm::NumericalProblem,
|
|
|
|
|
"Encountered a non-finite residual");
|
|
|
|
|
}
|
|
|
|
|
residualNorms.push_back(massBalanceResid);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
// the following residuals are not used in the oscillation detection now
|
|
|
|
|
const double wellFluxResid = detail::infinityNorm( residual_.well_flux_eq );
|
|
|
|
|
if (!std::isfinite(wellFluxResid)) {
|
|
|
|
|
OPM_THROW(Opm::NumericalProblem,
|
|
|
|
|
"Encountered a non-finite residual");
|
|
|
|
|
}
|
|
|
|
|
residualNorms.push_back(wellFluxResid);
|
|
|
|
|
|
|
|
|
|
const double wellResid = detail::infinityNorm( residual_.well_eq );
|
|
|
|
|
if (!std::isfinite(wellResid)) {
|
|
|
|
|
OPM_THROW(Opm::NumericalProblem,
|
|
|
|
|
"Encountered a non-finite residual");
|
|
|
|
|
}
|
|
|
|
|
residualNorms.push_back(wellResid);
|
|
|
|
|
|
|
|
|
|
return residualNorms;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
double
|
|
|
|
|
BlackoilPolymerModel<Grid>::convergenceReduction(const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases+1>& B,
|
2015-05-21 16:29:01 +02:00
|
|
|
const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases+1>& tempV,
|
|
|
|
|
const Eigen::Array<double, Eigen::Dynamic, MaxNumPhases+1>& R,
|
|
|
|
|
std::array<double,MaxNumPhases+1>& R_sum,
|
|
|
|
|
std::array<double,MaxNumPhases+1>& maxCoeff,
|
|
|
|
|
std::array<double,MaxNumPhases+1>& B_avg,
|
|
|
|
|
std::vector<double>& maxNormWell,
|
|
|
|
|
int nc,
|
|
|
|
|
int nw) const
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
|
|
|
|
// Do the global reductions
|
|
|
|
|
#if HAVE_MPI
|
|
|
|
|
if ( linsolver_.parallelInformation().type() == typeid(ParallelISTLInformation) )
|
|
|
|
|
{
|
|
|
|
|
const ParallelISTLInformation& info =
|
|
|
|
|
boost::any_cast<const ParallelISTLInformation&>(linsolver_.parallelInformation());
|
|
|
|
|
|
|
|
|
|
// Compute the global number of cells and porevolume
|
|
|
|
|
std::vector<int> v(nc, 1);
|
|
|
|
|
auto nc_and_pv = std::tuple<int, double>(0, 0.0);
|
|
|
|
|
auto nc_and_pv_operators = std::make_tuple(Opm::Reduction::makeGlobalSumFunctor<int>(),
|
|
|
|
|
Opm::Reduction::makeGlobalSumFunctor<double>());
|
|
|
|
|
auto nc_and_pv_containers = std::make_tuple(v, geo_.poreVolume());
|
|
|
|
|
info.computeReduction(nc_and_pv_containers, nc_and_pv_operators, nc_and_pv);
|
|
|
|
|
|
2015-05-25 00:06:17 +02:00
|
|
|
for ( int idx=0; idx<MaxNumPhases+1; ++idx )
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-25 00:10:11 +02:00
|
|
|
if (idx == MaxNumPhases || active_[idx]) { // Dealing with polymer *or* an active phase.
|
2015-05-19 14:40:38 +02:00
|
|
|
auto values = std::tuple<double,double,double>(0.0 ,0.0 ,0.0);
|
|
|
|
|
auto containers = std::make_tuple(B.col(idx),
|
|
|
|
|
tempV.col(idx),
|
|
|
|
|
R.col(idx));
|
|
|
|
|
auto operators = std::make_tuple(Opm::Reduction::makeGlobalSumFunctor<double>(),
|
|
|
|
|
Opm::Reduction::makeGlobalMaxFunctor<double>(),
|
|
|
|
|
Opm::Reduction::makeGlobalSumFunctor<double>());
|
|
|
|
|
info.computeReduction(containers, operators, values);
|
2015-05-21 16:29:01 +02:00
|
|
|
B_avg[idx] = std::get<0>(values)/std::get<0>(nc_and_pv);
|
|
|
|
|
maxCoeff[idx] = std::get<1>(values);
|
|
|
|
|
R_sum[idx] = std::get<2>(values);
|
2015-05-25 00:10:11 +02:00
|
|
|
if (idx != MaxNumPhases) { // We do not compute a well flux residual for polymer.
|
|
|
|
|
maxNormWell[idx] = 0.0;
|
|
|
|
|
for ( int w=0; w<nw; ++w ) {
|
|
|
|
|
maxNormWell[idx] = std::max(maxNormWell[idx], std::abs(residual_.well_flux_eq.value()[nw*idx + w]));
|
|
|
|
|
}
|
2015-05-21 16:29:01 +02:00
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
else
|
|
|
|
|
{
|
2015-05-21 16:29:01 +02:00
|
|
|
maxNormWell[idx] = R_sum[idx] = B_avg[idx] = maxCoeff[idx] = 0.0;
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
}
|
2015-05-25 00:06:17 +02:00
|
|
|
info.communicator().max(&maxNormWell[0], MaxNumPhases+1);
|
2015-05-19 14:40:38 +02:00
|
|
|
// Compute pore volume
|
2015-05-21 16:29:01 +02:00
|
|
|
#warning Missing polymer code for MPI version
|
2015-05-19 14:40:38 +02:00
|
|
|
return std::get<1>(nc_and_pv);
|
|
|
|
|
}
|
|
|
|
|
else
|
|
|
|
|
#endif
|
|
|
|
|
{
|
2015-05-25 00:06:17 +02:00
|
|
|
for ( int idx=0; idx<MaxNumPhases+1; ++idx )
|
2015-05-19 14:40:38 +02:00
|
|
|
{
|
2015-05-25 00:06:17 +02:00
|
|
|
if (idx == MaxNumPhases || active_[idx]) { // Dealing with polymer *or* an active phase.
|
2015-05-19 14:40:38 +02:00
|
|
|
B_avg[idx] = B.col(idx).sum()/nc;
|
|
|
|
|
maxCoeff[idx]=tempV.col(idx).maxCoeff();
|
|
|
|
|
R_sum[idx] = R.col(idx).sum();
|
|
|
|
|
}
|
|
|
|
|
else
|
|
|
|
|
{
|
|
|
|
|
R_sum[idx] = B_avg[idx] = maxCoeff[idx] =0.0;
|
|
|
|
|
}
|
2015-05-25 00:06:17 +02:00
|
|
|
if (idx != MaxNumPhases) { // We do not compute a well flux residual for polymer.
|
|
|
|
|
maxNormWell[idx] = 0.0;
|
|
|
|
|
for ( int w=0; w<nw; ++w ) {
|
|
|
|
|
maxNormWell[idx] = std::max(maxNormWell[idx], std::abs(residual_.well_flux_eq.value()[nw*idx + w]));
|
|
|
|
|
}
|
2015-05-21 16:29:01 +02:00
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
B_avg[MaxNumPhases] = B.col(MaxNumPhases).sum()/nc;
|
|
|
|
|
maxCoeff[MaxNumPhases]=tempV.col(MaxNumPhases).maxCoeff();
|
|
|
|
|
R_sum[MaxNumPhases] = R.col(MaxNumPhases).sum();
|
|
|
|
|
}
|
|
|
|
|
// Compute total pore volume
|
|
|
|
|
return geo_.poreVolume().sum();
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
|
2015-05-21 16:29:01 +02:00
|
|
|
|
|
|
|
|
|
|
|
|
|
|
2015-05-19 14:40:38 +02:00
|
|
|
template <class Grid>
|
|
|
|
|
bool
|
|
|
|
|
BlackoilPolymerModel<Grid>::getConvergence(const double dt, const int iteration)
|
|
|
|
|
{
|
|
|
|
|
const double tol_mb = param_.tolerance_mb_;
|
|
|
|
|
const double tol_cnv = param_.tolerance_cnv_;
|
|
|
|
|
const double tol_wells = param_.tolerance_wells_;
|
|
|
|
|
|
|
|
|
|
const int nc = Opm::AutoDiffGrid::numCells(grid_);
|
2015-05-21 16:29:01 +02:00
|
|
|
const int nw = wellsActive() ? wells().number_of_wells : 0;
|
2015-05-19 14:40:38 +02:00
|
|
|
const Opm::PhaseUsage& pu = fluid_.phaseUsage();
|
|
|
|
|
|
|
|
|
|
const V pv = geo_.poreVolume();
|
|
|
|
|
|
|
|
|
|
const std::vector<PhasePresence> cond = phaseCondition();
|
|
|
|
|
|
|
|
|
|
std::array<double,MaxNumPhases+1> CNV = {{0., 0., 0., 0.}};
|
|
|
|
|
std::array<double,MaxNumPhases+1> R_sum = {{0., 0., 0., 0.}};
|
|
|
|
|
std::array<double,MaxNumPhases+1> B_avg = {{0., 0., 0., 0.}};
|
|
|
|
|
std::array<double,MaxNumPhases+1> maxCoeff = {{0., 0., 0., 0.}};
|
|
|
|
|
std::array<double,MaxNumPhases+1> mass_balance_residual = {{0., 0., 0., 0.}};
|
2015-05-25 00:06:17 +02:00
|
|
|
std::array<double,MaxNumPhases> well_flux_residual = {{0., 0., 0.}};
|
2015-05-19 14:40:38 +02:00
|
|
|
std::size_t cols = MaxNumPhases+1; // needed to pass the correct type to Eigen
|
|
|
|
|
Eigen::Array<V::Scalar, Eigen::Dynamic, MaxNumPhases+1> B(nc, cols);
|
|
|
|
|
Eigen::Array<V::Scalar, Eigen::Dynamic, MaxNumPhases+1> R(nc, cols);
|
|
|
|
|
Eigen::Array<V::Scalar, Eigen::Dynamic, MaxNumPhases+1> tempV(nc, cols);
|
2015-05-21 16:29:01 +02:00
|
|
|
std::vector<double> maxNormWell(MaxNumPhases);
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
for ( int idx=0; idx<MaxNumPhases; ++idx )
|
|
|
|
|
{
|
|
|
|
|
if (active_[idx]) {
|
|
|
|
|
const int pos = pu.phase_pos[idx];
|
|
|
|
|
const ADB& tempB = rq_[pos].b;
|
|
|
|
|
B.col(idx) = 1./tempB.value();
|
|
|
|
|
R.col(idx) = residual_.material_balance_eq[idx].value();
|
|
|
|
|
tempV.col(idx) = R.col(idx).abs()/pv;
|
|
|
|
|
}
|
|
|
|
|
}
|
|
|
|
|
if (has_polymer_) {
|
|
|
|
|
const ADB& tempB = rq_[poly_pos_].b;
|
|
|
|
|
B.col(MaxNumPhases) = 1. / tempB.value();
|
|
|
|
|
R.col(MaxNumPhases) = residual_.material_balance_eq[poly_pos_].value();
|
|
|
|
|
tempV.col(MaxNumPhases) = R.col(MaxNumPhases).abs()/pv;
|
|
|
|
|
}
|
|
|
|
|
|
2015-05-21 16:29:01 +02:00
|
|
|
const double pvSum = convergenceReduction(B, tempV, R, R_sum, maxCoeff, B_avg,
|
|
|
|
|
maxNormWell, nc, nw);
|
2015-05-19 14:40:38 +02:00
|
|
|
|
|
|
|
|
bool converged_MB = true;
|
|
|
|
|
bool converged_CNV = true;
|
2015-05-21 16:29:01 +02:00
|
|
|
bool converged_Well = true;
|
2015-05-19 14:40:38 +02:00
|
|
|
// Finish computation
|
|
|
|
|
for ( int idx=0; idx<MaxNumPhases+1; ++idx )
|
|
|
|
|
{
|
|
|
|
|
CNV[idx] = B_avg[idx] * dt * maxCoeff[idx];
|
|
|
|
|
mass_balance_residual[idx] = std::abs(B_avg[idx]*R_sum[idx]) * dt / pvSum;
|
|
|
|
|
converged_MB = converged_MB && (mass_balance_residual[idx] < tol_mb);
|
|
|
|
|
converged_CNV = converged_CNV && (CNV[idx] < tol_cnv);
|
2015-05-25 00:06:17 +02:00
|
|
|
if (idx != MaxNumPhases) { // No well flux residual for polymer.
|
|
|
|
|
well_flux_residual[idx] = B_avg[idx] * dt * maxNormWell[idx];
|
|
|
|
|
converged_Well = converged_Well && (well_flux_residual[idx] < tol_wells);
|
|
|
|
|
}
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
2015-05-21 16:29:01 +02:00
|
|
|
const double residualWell = detail::infinityNormWell(residual_.well_eq,
|
|
|
|
|
linsolver_.parallelInformation());
|
|
|
|
|
converged_Well = converged_Well && (residualWell < Opm::unit::barsa);
|
2015-05-19 14:40:38 +02:00
|
|
|
const bool converged = converged_MB && converged_CNV && converged_Well;
|
|
|
|
|
|
|
|
|
|
// if one of the residuals is NaN, throw exception, so that the solver can be restarted
|
2015-05-25 00:06:17 +02:00
|
|
|
if (std::isnan(mass_balance_residual[Water]) || mass_balance_residual[Water] > maxResidualAllowed() ||
|
2015-05-19 14:40:38 +02:00
|
|
|
std::isnan(mass_balance_residual[Oil]) || mass_balance_residual[Oil] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(mass_balance_residual[Gas]) || mass_balance_residual[Gas] > maxResidualAllowed() ||
|
2015-05-25 00:06:17 +02:00
|
|
|
std::isnan(mass_balance_residual[MaxNumPhases]) || mass_balance_residual[MaxNumPhases] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(CNV[Water]) || CNV[Water] > maxResidualAllowed() ||
|
2015-05-19 14:40:38 +02:00
|
|
|
std::isnan(CNV[Oil]) || CNV[Oil] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(CNV[Gas]) || CNV[Gas] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(CNV[MaxNumPhases]) || CNV[MaxNumPhases] > maxResidualAllowed() ||
|
2015-05-21 16:29:01 +02:00
|
|
|
std::isnan(well_flux_residual[Water]) || well_flux_residual[Water] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(well_flux_residual[Oil]) || well_flux_residual[Oil] > maxResidualAllowed() ||
|
|
|
|
|
std::isnan(well_flux_residual[Gas]) || well_flux_residual[Gas] > maxResidualAllowed() ||
|
2015-05-19 14:40:38 +02:00
|
|
|
std::isnan(residualWell) || residualWell > maxResidualAllowed() )
|
|
|
|
|
{
|
2015-05-25 00:06:17 +02:00
|
|
|
OPM_THROW(Opm::NumericalProblem,"One of the residuals is NaN or too large!");
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
if ( terminal_output_ )
|
|
|
|
|
{
|
|
|
|
|
// Only rank 0 does print to std::cout
|
|
|
|
|
if (iteration == 0) {
|
2015-05-25 00:06:17 +02:00
|
|
|
std::cout << "\nIter MB(WATER) MB(OIL) MB(GAS) MB(POLY) CNVW CNVO CNVG CNVP W-FLUX(W) W-FLUX(O) W-FLUX(G)\n";
|
2015-05-19 14:40:38 +02:00
|
|
|
}
|
|
|
|
|
const std::streamsize oprec = std::cout.precision(3);
|
|
|
|
|
const std::ios::fmtflags oflags = std::cout.setf(std::ios::scientific);
|
|
|
|
|
std::cout << std::setw(4) << iteration
|
|
|
|
|
<< std::setw(11) << mass_balance_residual[Water]
|
|
|
|
|
<< std::setw(11) << mass_balance_residual[Oil]
|
|
|
|
|
<< std::setw(11) << mass_balance_residual[Gas]
|
|
|
|
|
<< std::setw(11) << mass_balance_residual[MaxNumPhases]
|
|
|
|
|
<< std::setw(11) << CNV[Water]
|
|
|
|
|
<< std::setw(11) << CNV[Oil]
|
|
|
|
|
<< std::setw(11) << CNV[Gas]
|
|
|
|
|
<< std::setw(11) << CNV[MaxNumPhases]
|
2015-05-21 16:29:01 +02:00
|
|
|
<< std::setw(11) << well_flux_residual[Water]
|
|
|
|
|
<< std::setw(11) << well_flux_residual[Oil]
|
|
|
|
|
<< std::setw(11) << well_flux_residual[Gas]
|
2015-05-19 14:40:38 +02:00
|
|
|
<< std::endl;
|
|
|
|
|
std::cout.precision(oprec);
|
|
|
|
|
std::cout.flags(oflags);
|
|
|
|
|
}
|
|
|
|
|
return converged;
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
template <class Grid>
|
|
|
|
|
ADB
|
|
|
|
|
BlackoilPolymerModel<Grid>::computeMc(const SolutionState& state) const
|
|
|
|
|
{
|
|
|
|
|
return polymer_props_ad_.polymerWaterVelocityRatio(state.concentration);
|
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
} // namespace Opm
|
|
|
|
|
|
|
|
|
|
#endif // OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
|