fix some bugs.

This commit is contained in:
Liu Ming 2014-10-09 10:04:17 +08:00
parent be927741b8
commit 73c2032974
4 changed files with 50 additions and 59 deletions

View File

@ -197,8 +197,9 @@ try
simtimer.init(timeMap);
//Check for WPOLYMER presence in last report step to decide
//polymer injection control type.
const bool polymer = deck->hasKeywrod("POLYMER");
const bool polymer = deck->hasKeyword("POLYMER");
const bool use_wpolymer = deck->hasKeyword("WPOLYMER");
std::shared_ptr<PolymerInflowInterface> polymer_inflow;
if (polymer){
if (!use_wpolymer) {
OPM_MESSAGE("Warning: simulate polymer injection without WPOLYMER.");
@ -207,7 +208,6 @@ try
OPM_MESSAGE("Warning: Using WPOLYMER to control injection since it was found in deck."
"You seem to be trying to control it via parameter poly_start_days (etc.) as well.");
}
std::shared_ptr<PolymerInflowInterface> polymer_inflow;
}
} else {
if (use_wpolymer) {

View File

@ -68,17 +68,17 @@ namespace Opm {
/// \param[in] rock_comp_props if non-null, rock compressibility properties
/// \param[in] wells well structure
/// \param[in] linsolver linear solver
FullyImplicitBlackoilSolver(const parameter::ParameterGroup& 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);
FullyImplicitBlackoilPolymerSolver(const parameter::ParameterGroup& 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);
/// \brief Set threshold pressures that prevent or reduce flow.
/// This prevents flow across faces if the potential
@ -251,18 +251,16 @@ namespace Opm {
void
computeMassFlux(const V& trans,
const ADB& mc,
const ADB& kro,
const ADB& krw_eff,
const ADB& krg,
const std::vector<ADB>& kr,
const std::vector<ADB>& phasePressure,
const SolutionState& state);
void
computeCmax(PolymerBlackoilState& state,
const ADB& c);
computeCmax(PolymerBlackoilState& state,
const ADB& c);
ADB
computeMc(const SolutionState& state) const;
computeMc(const SolutionState& state) const;
ADB
rockPorosity(const ADB& p) const;
@ -320,8 +318,6 @@ namespace Opm {
const ADB& so,
const std::vector<int>& cells) const;
ADB
computeMc(const SolutionState& state) const;
ADB
poroMult(const ADB& p) const;
@ -371,6 +367,6 @@ namespace Opm {
};
} // namespace Opm
#include "FullyImplicitBlackoilSolver_impl.hpp"
#include "FullyImplicitBlackoilPolymerSolver_impl.hpp"
#endif // OPM_FULLYIMPLICITBLACKOILSOLVER_HEADER_INCLUDED
#endif // OPM_FULLYIMPLICITBLACKOILPOLYMERSOLVER_HEADER_INCLUDED

View File

@ -264,7 +264,8 @@ namespace {
FullyImplicitBlackoilPolymerSolver<T>::
step(const double dt,
PolymerBlackoilState& x ,
WellStateFullyImplicitBlackoil& xw)
WellStateFullyImplicitBlackoil& xw,
const std::vector<double>& polymer_inflow)
{
const V pvdt = geo_.poreVolume() / dt;
@ -625,7 +626,6 @@ namespace {
// 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;
}
}
@ -633,7 +633,7 @@ namespace {
template<class T>
void FullyImplicitBlackoilPolymerSovler<T>::computeCmax(PolymerBlackoilState& state,
void FullyImplicitBlackoilPolymerSolver<T>::computeCmax(PolymerBlackoilState& state,
const ADB& c)
{
const int nc = numCells(grid_);
@ -724,7 +724,7 @@ namespace {
void
FullyImplicitBlackoilPolymerSolver<T>::
assemble(const V& pvdt,
const BlackoilState& x ,
const PolymerBlackoilState& x ,
WellStateFullyImplicitBlackoil& xw,
const std::vector<double>& polymer_inflow)
{
@ -798,7 +798,7 @@ namespace {
const Opm::PhaseUsage& pu = fluid_.phaseUsage();
// Add polymer equation.
if (has_polymer_) {
residula_.material_balance_eq[poly_pos_] = pvdt*(rq_[poly_pos_].accum[1] - rq_[poly_pos_].accum[0])
residual_.material_balance_eq[poly_pos_] = pvdt*(rq_[poly_pos_].accum[1] - rq_[poly_pos_].accum[0])
+ ops_.div*rq_[poly_pos_].mflux;
}
@ -1001,7 +1001,7 @@ namespace {
}
// Add well contributions to polymer mass_balance equation
if (has_polmer_) {
if (has_polymer_) {
const ADB mc = computeMc(state);
const V polyin = Eigen::Map<const V>(&polymer_inflow[0], nc);
const V poly_in_perf = subset(polyin, well_cells);
@ -1533,6 +1533,8 @@ namespace {
//Polymer concentration updates.
if (has_polymer_) {
const V dc = subset(dx, Span(nc));
int varstart = nc;
const V c_old = Eigen::Map<const V>(&state.concentration()[0], nc, 1);
const V c = (c_old - dc).max(zero);
std::copy(&c[0], &c[0] + nc, state.concentration().begin());
@ -1679,13 +1681,6 @@ namespace {
const std::vector<ADB>& phasePressure,
const SolutionState& state)
{
if (has_polymer_) {
const ADB cmax = ADB::constant(cmax_, state.concentration.blockPattern());
ADB krw_eff = ADB::null();
const ADB mc = computeMc(state);
ADB inv_wat_eff_vis = ADB::null();
}
for (int phase = 0; phase < fluid_.numPhases(); ++phase) {
const int canonicalPhaseIdx = canph_[phase];
const std::vector<PhasePresence> cond = phaseCondition();
@ -1693,14 +1688,21 @@ namespace {
const ADB tr_mult = transMult(state.pressure);
const ADB mu = fluidViscosity(canonicalPhaseIdx, phasePressure[canonicalPhaseIdx], state.rs, state.rv, cond, cells_);
rq_[canonicalPhaseIdx].mob = tr_mult * kr[canonicalPhaseIdx] / mu;
if (cononicalPhaseIdx == Water) {
if (canonicalPhaseIdx == Water) {
if(has_polymer_) {
krw_eff = polymer_props_ad_.effectiveRelPerm(state.concentration,
cmax,
kr[cononicalPhaseIdx],
state.saturation[cononicalPhaseIdx]);
inv_wat_eff_vis = polymer_props_ad_.effectiveInvWaterVisc(state.concentration, mu.value().data());
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[canonicalPhaseIdx],
state.saturation[canonicalPhaseIdx]);
ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(state.concentration, mu.value().data());
rq_[canonicalPhaseIdx].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_[canph_[Water]].b;
rq_[poly_pos_].head = rq_[canph_[Water]].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;
}
}
@ -1726,14 +1728,6 @@ namespace {
const ADB& mob = rq_[canonicalPhaseIdx].mob;
rq_[canonicalPhaseIdx].mflux = upwind.select(b * mob) * head;
}
if (has_polymer_) {
rq_[poly_pos_].mob = tr_mult * mc * krw_eff * inv_wat_eff_visc;
rq_[poly_pos_].b = rq_[canph_[Water]].b;
rq_[poly_pos_].head = rq_[canph_[Water]].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;
}
}
@ -2119,7 +2113,7 @@ namespace {
template<class T>
ADB
FullyImplicitBlackoilPolymerSovler<T>::computeMc(const SolutionState& state) const
FullyImplicitBlackoilPolymerSolver<T>::computeMc(const SolutionState& state) const
{
return polymer_props_ad_.polymerWaterVelocityRatio(state.concentration);
}
@ -2146,7 +2140,7 @@ namespace {
jacs[block] = dpm_diag * p.derivative()[block];
}
return ADB::function(pm, jacs);
else {
} else {
return ADB::constant(V::Constant(n, 1.0), p.blockPattern());
}
}

View File

@ -20,8 +20,9 @@
#include <opm/polymer/fullyimplicit/SimulatorFullyImplicitBlackoilPolymerOutput.hpp>
#include <opm/polymer/fullyimplicit/SimulatorFullyImplicitBlackoilPolymer.hpp>
#include <opm/polymer/fullyimplicit/FullyImplicitBlackoilSolverPolymer.hpp>
#include <opm/polymer/fullyimplicit/FullyImplicitBlackoilPolymerSolver.hpp>
#include <opm/polymer/PolymerBlackoilState.hpp>
#include <opm/polymer/PolymerInflow.hpp>
#include <opm/core/utility/parameters/ParameterGroup.hpp>
#include <opm/core/utility/ErrorMacros.hpp>
@ -115,7 +116,7 @@ namespace Opm
// Solvers
const DerivedGeology& geo_;
NewtonIterationBlackoilInterface& solver_;
std::shared_ptr<PolymerInflowInterface> polymer_inflow,
std::shared_ptr<PolymerInflowInterface> polymer_inflow_;
// Misc. data
std::vector<int> allcells_;
const bool has_disgas_;
@ -237,7 +238,7 @@ namespace Opm
const Grid& grid,
const DerivedGeology& geo,
BlackoilPropsAdInterface& props,
const PolymerPropdAd& polymer_props,
const PolymerPropsAd& polymer_props,
const RockCompressibility* rock_comp_props,
NewtonIterationBlackoilInterface& linsolver,
std::shared_ptr<PolymerInflowInterface> polymer_inflow,
@ -355,10 +356,10 @@ namespace Opm
computeRESV(timer.currentStepNum(), wells, state, well_state);
// compute polymer inflow
std::vector<double> polymer_inflow_c(Opm::UgGridHelpers::numCells(grid_));
if (has_polymer_) {
std::vector<double> polymer_inflow_c(Opm::UgGridHelpers::numCells(grid_));
polymer_inflow_.reset(new PolymerFromDeck(deck, wells, Opm::UgGridHelpers::numCells(grid_)));
polymer_inflow_.getInflowValues(timer.simulationTimeElapsed(),
polymer_inflow_.reset(new PolymerInflowFromDeck(deck, wells, Opm::UgGridHelpers::numCells(grid_)));
*polymer_inflow_.getInflowValues(timer.simulationTimeElapsed(),
timer.simulationTimeElapsed() + timer.currentStepLength(),
polymer_inflow_c);
}