mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
putting update of water injection mobility to a separate function
computeInjectionMobility() in BlackoilPolymerModel
This commit is contained in:
parent
4f3ce37b8e
commit
b8e0cda83a
@ -227,6 +227,10 @@ namespace Opm {
|
|||||||
computeAccum(const SolutionState& state,
|
computeAccum(const SolutionState& state,
|
||||||
const int aix );
|
const int aix );
|
||||||
|
|
||||||
|
void
|
||||||
|
computeInjectionMobility(const SolutionState& state,
|
||||||
|
std::vector<ADB>& mob_perfcells);
|
||||||
|
|
||||||
void
|
void
|
||||||
assembleMassBalanceEq(const SolutionState& state);
|
assembleMassBalanceEq(const SolutionState& state);
|
||||||
|
|
||||||
|
@ -540,58 +540,9 @@ namespace Opm {
|
|||||||
std::vector<ADB> b_perfcells;
|
std::vector<ADB> b_perfcells;
|
||||||
wellModel().extractWellPerfProperties(state, sd_.rq, mob_perfcells, b_perfcells);
|
wellModel().extractWellPerfProperties(state, sd_.rq, mob_perfcells, b_perfcells);
|
||||||
|
|
||||||
// Calculating the mobility for the polymer injection peforations
|
// updating the the injection mobility related to polymer injection when necessary
|
||||||
if (has_polymer_ && wellModel().localWellsActive()) {
|
// only the mobility of water phase is updated
|
||||||
const std::vector<int> well_cells = wellModel().wellOps().well_cells;
|
computeInjectionMobility(state, mob_perfcells);
|
||||||
const int nperf = well_cells.size();
|
|
||||||
|
|
||||||
// Calculating the drawdown to decide the injection perforation
|
|
||||||
const ADB& p_perfcells = subset(state.pressure, well_cells);
|
|
||||||
const V& cdp = wellModel().wellPerforationPressureDiffs();
|
|
||||||
const ADB perfpressure = (wellModel().wellOps().w2p * state.bhp) + cdp;
|
|
||||||
// Pressure drawdown (also used to determine direction of flow)
|
|
||||||
const ADB drawdown = p_perfcells - perfpressure;
|
|
||||||
|
|
||||||
// Polymer concentration in the perforations
|
|
||||||
const ADB c_perfcells = subset(state.concentration, well_cells);
|
|
||||||
|
|
||||||
// Distinguishing the injection perforation from other perforation
|
|
||||||
// The value is the location in the well_cell array, not the global index
|
|
||||||
std::vector<int> polymer_inj_cells;
|
|
||||||
std::vector<int> other_well_cells;
|
|
||||||
|
|
||||||
polymer_inj_cells.reserve(nperf);
|
|
||||||
other_well_cells.reserve(nperf);
|
|
||||||
|
|
||||||
for (int c = 0; c < nperf; ++c) {
|
|
||||||
// TODO: more tests need to be done for this criterion
|
|
||||||
if (drawdown.value()[c] < 0.0 && c_perfcells.value()[c] > 0.0) {
|
|
||||||
polymer_inj_cells.push_back(c);
|
|
||||||
} else {
|
|
||||||
other_well_cells.push_back(c);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
|
|
||||||
// there is some polymer injection process going
|
|
||||||
if ( !polymer_inj_cells.empty() ) {
|
|
||||||
// the mobility need to be recalculated for the polymer injection cells
|
|
||||||
const int water_pos = fluid_.phaseUsage().phase_pos[Water];
|
|
||||||
const ADB mu_perfcells = subset(sd_.rq[water_pos].mu, well_cells);
|
|
||||||
|
|
||||||
const ADB c_poly_inj_cells = subset(c_perfcells, polymer_inj_cells);
|
|
||||||
const ADB mu_poly_inj_cells = subset(mu_perfcells, polymer_inj_cells); // water viscosity
|
|
||||||
|
|
||||||
const ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(c_poly_inj_cells, mu_poly_inj_cells.value());
|
|
||||||
const ADB fully_mixing_visc = polymer_props_ad_.viscMult(c_poly_inj_cells) * mu_poly_inj_cells;
|
|
||||||
|
|
||||||
// the original mobility for the polymer injection well cells
|
|
||||||
ADB mob_polymer_inj = subset(mob_perfcells[water_pos], polymer_inj_cells);
|
|
||||||
const ADB mob_others = subset(mob_perfcells[water_pos], other_well_cells);
|
|
||||||
|
|
||||||
mob_polymer_inj = mob_polymer_inj / inv_wat_eff_visc / fully_mixing_visc;
|
|
||||||
mob_perfcells[water_pos] = superset(mob_polymer_inj, polymer_inj_cells, nperf) + superset(mob_others, other_well_cells, nperf);
|
|
||||||
}
|
|
||||||
}
|
|
||||||
|
|
||||||
if (param_.solve_welleq_initially_ && initial_assembly) {
|
if (param_.solve_welleq_initially_ && initial_assembly) {
|
||||||
// solve the well equations as a pre-processing step
|
// solve the well equations as a pre-processing step
|
||||||
@ -838,6 +789,71 @@ namespace Opm {
|
|||||||
return;
|
return;
|
||||||
}
|
}
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
|
||||||
|
template<class Grid>
|
||||||
|
void
|
||||||
|
BlackoilPolymerModel<Grid>::
|
||||||
|
computeInjectionMobility(const SolutionState& state,
|
||||||
|
std::vector<ADB>& mob_perfcells)
|
||||||
|
{
|
||||||
|
// Calculating the mobility for the polymer injection peforations
|
||||||
|
if (has_polymer_ && wellModel().localWellsActive()) {
|
||||||
|
const std::vector<int> well_cells = wellModel().wellOps().well_cells;
|
||||||
|
const int nperf = well_cells.size();
|
||||||
|
|
||||||
|
// Calculating the drawdown to decide the injection perforation
|
||||||
|
const ADB& p_perfcells = subset(state.pressure, well_cells);
|
||||||
|
const V& cdp = wellModel().wellPerforationPressureDiffs();
|
||||||
|
const ADB perfpressure = (wellModel().wellOps().w2p * state.bhp) + cdp;
|
||||||
|
// Pressure drawdown (also used to determine direction of flow)
|
||||||
|
const ADB drawdown = p_perfcells - perfpressure;
|
||||||
|
|
||||||
|
// Polymer concentration in the perforations
|
||||||
|
const ADB c_perfcells = subset(state.concentration, well_cells);
|
||||||
|
|
||||||
|
// Distinguishing the injection perforation from other perforation
|
||||||
|
// The value is the location in the well_cell array, not the global index
|
||||||
|
std::vector<int> polymer_inj_cells;
|
||||||
|
std::vector<int> other_well_cells;
|
||||||
|
|
||||||
|
polymer_inj_cells.reserve(nperf);
|
||||||
|
other_well_cells.reserve(nperf);
|
||||||
|
|
||||||
|
for (int c = 0; c < nperf; ++c) {
|
||||||
|
// TODO: more tests need to be done for this criterion
|
||||||
|
if (drawdown.value()[c] < 0.0 && c_perfcells.value()[c] > 0.0) {
|
||||||
|
polymer_inj_cells.push_back(c);
|
||||||
|
} else {
|
||||||
|
other_well_cells.push_back(c);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
// there is some polymer injection process going
|
||||||
|
if ( !polymer_inj_cells.empty() ) {
|
||||||
|
// the mobility need to be recalculated for the polymer injection cells
|
||||||
|
const int water_pos = fluid_.phaseUsage().phase_pos[Water];
|
||||||
|
const ADB mu_perfcells = subset(sd_.rq[water_pos].mu, well_cells);
|
||||||
|
|
||||||
|
const ADB c_poly_inj_cells = subset(c_perfcells, polymer_inj_cells);
|
||||||
|
const ADB mu_poly_inj_cells = subset(mu_perfcells, polymer_inj_cells); // water viscosity
|
||||||
|
|
||||||
|
const ADB inv_wat_eff_visc = polymer_props_ad_.effectiveInvWaterVisc(c_poly_inj_cells, mu_poly_inj_cells.value());
|
||||||
|
const ADB fully_mixing_visc = polymer_props_ad_.viscMult(c_poly_inj_cells) * mu_poly_inj_cells;
|
||||||
|
|
||||||
|
// the original mobility for the polymer injection well cells
|
||||||
|
ADB mob_polymer_inj = subset(mob_perfcells[water_pos], polymer_inj_cells);
|
||||||
|
const ADB mob_others = subset(mob_perfcells[water_pos], other_well_cells);
|
||||||
|
|
||||||
|
mob_polymer_inj = mob_polymer_inj / inv_wat_eff_visc / fully_mixing_visc;
|
||||||
|
mob_perfcells[water_pos] = superset(mob_polymer_inj, polymer_inj_cells, nperf) +
|
||||||
|
superset(mob_others, other_well_cells, nperf);
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
} // namespace Opm
|
} // namespace Opm
|
||||||
|
|
||||||
#endif // OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
|
#endif // OPM_BLACKOILPOLYMERMODEL_IMPL_HEADER_INCLUDED
|
||||||
|
Loading…
Reference in New Issue
Block a user