mirror of
https://github.com/OPM/opm-simulators.git
synced 2025-02-25 18:55:30 -06:00
Add saturation computation to and rename computer class.
Opm::equil::DeckDependent::PhasePressureComputer -> Opm::equil::DeckDependent::PhasePressureSaturationComputer
This commit is contained in:
parent
f8f9dc538a
commit
0a8ae66046
@ -871,42 +871,47 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
|
|
||||||
template <class InputDeck>
|
template <class InputDeck>
|
||||||
class PhasePressureComputer;
|
class PhasePressureSaturationComputer;
|
||||||
|
|
||||||
template <>
|
template <>
|
||||||
class PhasePressureComputer<Opm::EclipseGridParser> {
|
class PhasePressureSaturationComputer<Opm::EclipseGridParser> {
|
||||||
public:
|
public:
|
||||||
PhasePressureComputer(const BlackoilPropertiesInterface& props,
|
PhasePressureSaturationComputer(const BlackoilPropertiesInterface& props,
|
||||||
const EclipseGridParser& deck ,
|
const EclipseGridParser& deck ,
|
||||||
const UnstructuredGrid& G ,
|
const UnstructuredGrid& G ,
|
||||||
const double grav = unit::gravity)
|
const double grav = unit::gravity)
|
||||||
: pp_(props.numPhases(),
|
: pp_(props.numPhases(),
|
||||||
|
std::vector<double>(G.number_of_cells)),
|
||||||
|
sat_(props.numPhases(),
|
||||||
std::vector<double>(G.number_of_cells))
|
std::vector<double>(G.number_of_cells))
|
||||||
{
|
{
|
||||||
const std::vector<EquilRecord> rec = getEquil(deck);
|
const std::vector<EquilRecord> rec = getEquil(deck);
|
||||||
const RegionMapping<> eqlmap(equilnum(deck, G));
|
const RegionMapping<> eqlmap(equilnum(deck, G));
|
||||||
|
|
||||||
calcII(eqlmap, rec, props, G, grav);
|
calcPressII(eqlmap, rec, props, G, grav);
|
||||||
|
calcSat(eqlmap, rec, props, G, grav);
|
||||||
}
|
}
|
||||||
|
|
||||||
typedef std::vector<double> PVal;
|
typedef std::vector<double> PVal;
|
||||||
typedef std::vector<PVal> PPress;
|
typedef std::vector<PVal> PPress;
|
||||||
|
|
||||||
const PPress& press() const { return pp_; }
|
const PPress& press() const { return pp_; }
|
||||||
|
const PPress& saturation() const { return sat_; }
|
||||||
|
|
||||||
private:
|
private:
|
||||||
typedef DensityCalculator<BlackoilPropertiesInterface> RhoCalc;
|
typedef DensityCalculator<BlackoilPropertiesInterface> RhoCalc;
|
||||||
typedef EquilReg<RhoCalc> EqReg;
|
typedef EquilReg<RhoCalc> EqReg;
|
||||||
|
|
||||||
PPress pp_;
|
PPress pp_;
|
||||||
|
PPress sat_;
|
||||||
|
|
||||||
template <class RMap>
|
template <class RMap>
|
||||||
void
|
void
|
||||||
calcII(const RMap& reg ,
|
calcPressII(const RMap& reg ,
|
||||||
const std::vector< EquilRecord >& rec ,
|
const std::vector< EquilRecord >& rec ,
|
||||||
const Opm::BlackoilPropertiesInterface& props,
|
const Opm::BlackoilPropertiesInterface& props,
|
||||||
const UnstructuredGrid& G ,
|
const UnstructuredGrid& G ,
|
||||||
const double grav)
|
const double grav)
|
||||||
{
|
{
|
||||||
typedef miscibility::NoMixing NoMix;
|
typedef miscibility::NoMixing NoMix;
|
||||||
|
|
||||||
@ -937,6 +942,58 @@ namespace Opm
|
|||||||
}
|
}
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
|
template <class RMap>
|
||||||
|
void
|
||||||
|
calcSat(const RMap& reg ,
|
||||||
|
const std::vector< EquilRecord >& rec ,
|
||||||
|
const Opm::BlackoilPropertiesInterface& props,
|
||||||
|
const UnstructuredGrid& G ,
|
||||||
|
const double grav)
|
||||||
|
{
|
||||||
|
typedef miscibility::NoMixing NoMix;
|
||||||
|
|
||||||
|
for (typename RMap::RegionId
|
||||||
|
r = 0, nr = reg.numRegions();
|
||||||
|
r < nr; ++r)
|
||||||
|
{
|
||||||
|
const typename RMap::CellRange cells = reg.cells(r);
|
||||||
|
|
||||||
|
const int repcell = *cells.begin();
|
||||||
|
const RhoCalc calc(props, repcell);
|
||||||
|
|
||||||
|
const EqReg eqreg(rec[r], calc, NoMix(), NoMix(),
|
||||||
|
props.phaseUsage());
|
||||||
|
|
||||||
|
const PPress press = phasePressures(G, eqreg, cells, grav);
|
||||||
|
const PPress sat = phaseSaturations(eqreg, cells, props, press);
|
||||||
|
|
||||||
|
for (int p = 0, np = props.numPhases(); p < np; ++p) {
|
||||||
|
PVal& d = pp_[p];
|
||||||
|
PVal::const_iterator s = press[p].begin();
|
||||||
|
for (typename RMap::CellRange::const_iterator
|
||||||
|
c = cells.begin(),
|
||||||
|
e = cells.end();
|
||||||
|
c != e; ++c, ++s)
|
||||||
|
{
|
||||||
|
d[*c] = *s;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
for (int p = 0, np = props.numPhases(); p < np; ++p) {
|
||||||
|
PVal& d = sat_[p];
|
||||||
|
PVal::const_iterator s = sat[p].begin();
|
||||||
|
for (typename RMap::CellRange::const_iterator
|
||||||
|
c = cells.begin(),
|
||||||
|
e = cells.end();
|
||||||
|
c != e; ++c, ++s)
|
||||||
|
{
|
||||||
|
d[*c] = *s;
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
}
|
||||||
|
|
||||||
|
|
||||||
};
|
};
|
||||||
} // namespace DeckDependent
|
} // namespace DeckDependent
|
||||||
} // namespace equil
|
} // namespace equil
|
||||||
|
Loading…
Reference in New Issue
Block a user