22 #ifndef OPM_INITSTATEEQUIL_HEADER_INCLUDED
23 #define OPM_INITSTATEEQUIL_HEADER_INCLUDED
32 #include <opm/parser/eclipse/Utility/EquilWrapper.hpp>
33 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
67 const BlackoilPropertiesInterface& props,
68 const Opm::DeckConstPtr deck,
69 const Opm::EclipseStateConstPtr eclipseState,
71 BlackoilState& state);
121 template <
class Gr
id,
class Region,
class CellRange>
122 std::vector< std::vector<double> >
125 const CellRange& cells,
154 template <
class Gr
id,
class Region,
class CellRange>
155 std::vector< std::vector<double> >
158 const CellRange& cells,
159 BlackoilPropertiesInterface& props,
160 const std::vector<double> swat_init,
161 std::vector< std::vector<double> >& phase_pressures);
183 template <
class Gr
id,
class CellRangeType>
185 const CellRangeType& cells,
186 const std::vector<double> oil_pressure,
188 const Miscibility::RsFunction& rs_func,
189 const std::vector<double> gas_saturation);
191 namespace DeckDependent {
193 std::vector<EquilRecord>
194 getEquil(
const Opm::DeckConstPtr deck)
196 if (deck->hasKeyword(
"EQUIL")) {
198 Opm::EquilWrapper eql(deck->getKeyword(
"EQUIL"));
200 const int nrec = eql.numRegions();
202 std::vector<EquilRecord> ret;
204 for (
int r = 0; r < nrec; ++r) {
208 { eql.datumDepth(r) ,
209 eql.datumDepthPressure(r) }
211 { eql.waterOilContactDepth(r) ,
212 eql.waterOilContactCapillaryPressure(r) }
214 { eql.gasOilContactDepth(r) ,
215 eql.gasOilContactCapillaryPressure(r) }
217 eql.liveOilInitProceedure(r)
219 eql.wetGasInitProceedure(r)
221 eql.initializationTargetAccuracy(r)
224 OPM_THROW(std::domain_error,
225 "kw EQUIL, item 9: Only N=0 supported.");
227 ret.push_back(record);
233 OPM_THROW(std::domain_error,
234 "Deck does not provide equilibration data.");
241 equilnum(
const Opm::DeckConstPtr deck,
242 const Opm::EclipseStateConstPtr eclipseState,
245 std::vector<int> eqlnum;
246 if (deck->hasKeyword(
"EQLNUM")) {
249 const std::vector<int>& e =
250 eclipseState->getIntGridProperty(
"EQLNUM")->getData();
252 for (
int cell = 0; cell < nc; ++cell) {
253 const int deck_pos = (gc == NULL) ? cell : gc[cell];
254 eqlnum[cell] = e[deck_pos] - 1;
267 class InitialStateComputer {
270 InitialStateComputer(BlackoilPropertiesInterface& props,
271 const Opm::DeckConstPtr deck,
272 const Opm::EclipseStateConstPtr eclipseState,
275 : pp_(props.numPhases(),
277 sat_(props.numPhases(),
283 const std::vector<EquilRecord> rec = getEquil(deck);
284 std::shared_ptr<const TableManager> tables = eclipseState->getTableManager();
286 const RegionMapping<> eqlmap(equilnum(deck, eclipseState, G));
289 rs_func_.reserve(rec.size());
290 if (deck->hasKeyword(
"DISGAS")) {
291 const TableContainer& rsvdTables = tables->getRsvdTables();
292 for (
size_t i = 0; i < rec.size(); ++i) {
293 const int cell = *(eqlmap.cells(i).begin());
294 if (rec[i].live_oil_table_index > 0) {
295 if (rsvdTables.size() > 0 && size_t(rec[i].live_oil_table_index) <= rsvdTables.size()) {
296 const RsvdTable& rsvdTable = rsvdTables.getTable<RsvdTable>(i);
297 rs_func_.push_back(std::make_shared<Miscibility::RsVD>(props,
299 rsvdTable.getDepthColumn(),
300 rsvdTable.getRsColumn()));
302 OPM_THROW(std::runtime_error,
"Cannot initialise: RSVD table " << (rec[i].live_oil_table_index) <<
" not available.");
305 if (rec[i].goc.depth != rec[i].main.depth) {
306 OPM_THROW(std::runtime_error,
307 "Cannot initialise: when no explicit RSVD table is given, \n"
308 "datum depth must be at the gas-oil-contact. "
309 "In EQUIL region " << (i + 1) <<
" (counting from 1), this does not hold.");
311 const double p_contact = rec[i].main.press;
312 const double T_contact = 273.15 + 20;
313 rs_func_.push_back(std::make_shared<Miscibility::RsSatAtContact>(props, cell, p_contact, T_contact));
317 for (
size_t i = 0; i < rec.size(); ++i) {
318 rs_func_.push_back(std::make_shared<Miscibility::NoMixing>());
322 rv_func_.reserve(rec.size());
323 if (deck->hasKeyword(
"VAPOIL")) {
324 const TableContainer& rvvdTables = tables->getRvvdTables();
325 for (
size_t i = 0; i < rec.size(); ++i) {
326 const int cell = *(eqlmap.cells(i).begin());
327 if (rec[i].wet_gas_table_index > 0) {
328 if (rvvdTables.size() > 0 && size_t(rec[i].wet_gas_table_index) <= rvvdTables.size()) {
329 const RvvdTable& rvvdTable = rvvdTables.getTable<RvvdTable>(i);
330 rv_func_.push_back(std::make_shared<Miscibility::RvVD>(props,
332 rvvdTable.getDepthColumn(),
333 rvvdTable.getRvColumn()));
335 OPM_THROW(std::runtime_error,
"Cannot initialise: RVVD table " << (rec[i].wet_gas_table_index) <<
" not available.");
338 if (rec[i].goc.depth != rec[i].main.depth) {
339 OPM_THROW(std::runtime_error,
340 "Cannot initialise: when no explicit RVVD table is given, \n"
341 "datum depth must be at the gas-oil-contact. "
342 "In EQUIL region " << (i + 1) <<
" (counting from 1), this does not hold.");
344 const double p_contact = rec[i].main.press + rec[i].goc.press;
345 const double T_contact = 273.15 + 20;
346 rv_func_.push_back(std::make_shared<Miscibility::RvSatAtContact>(props, cell, p_contact, T_contact));
350 for (
size_t i = 0; i < rec.size(); ++i) {
351 rv_func_.push_back(std::make_shared<Miscibility::NoMixing>());
357 if (deck->hasKeyword(
"SWATINIT")) {
358 const std::vector<double>& swat_init = eclipseState->getDoubleGridProperty(
"SWATINIT")->getData();
360 swat_init_.resize(nc);
362 for (
int c = 0; c < nc; ++c) {
363 const int deck_pos = (gc == NULL) ? c : gc[c];
364 swat_init_[c] = swat_init[deck_pos];
369 calcPressSatRsRv(eqlmap, rec, props, G, grav);
375 typedef std::vector<double> Vec;
376 typedef std::vector<Vec> PVec;
378 const PVec& press()
const {
return pp_; }
379 const PVec& saturation()
const {
return sat_; }
380 const Vec& rs()
const {
return rs_; }
381 const Vec& rv()
const {
return rv_; }
384 typedef DensityCalculator<BlackoilPropertiesInterface> RhoCalc;
385 typedef EquilReg<RhoCalc> EqReg;
387 std::vector< std::shared_ptr<Miscibility::RsFunction> > rs_func_;
388 std::vector< std::shared_ptr<Miscibility::RsFunction> > rv_func_;
396 template <
class RMap,
class Gr
id>
398 calcPressSatRsRv(
const RMap& reg ,
399 const std::vector< EquilRecord >& rec ,
404 for (
const auto& r : reg.activeRegions()) {
405 const auto& cells = reg.cells(r);
406 const int repcell = *cells.begin();
408 const RhoCalc calc(props, repcell);
409 const EqReg eqreg(rec[r], calc,
410 rs_func_[r], rv_func_[r],
414 const std::vector<double>& temp =
temperature(G, eqreg, cells);
416 const PVec sat =
phaseSaturations(G, eqreg, cells, props, swat_init_, pressures);
419 for (
int p = 0; p < np; ++p) {
420 copyFromRegion(pressures[p], cells, pp_[p]);
421 copyFromRegion(sat[p], cells, sat_[p]);
427 const Vec rs_vals =
computeRs(G, cells, pressures[oilpos], temp, *(rs_func_[r]), sat[gaspos]);
428 const Vec rv_vals =
computeRs(G, cells, pressures[gaspos], temp, *(rv_func_[r]), sat[oilpos]);
429 copyFromRegion(rs_vals, cells, rs_);
430 copyFromRegion(rv_vals, cells, rv_);
435 template <
class CellRangeType>
436 void copyFromRegion(
const Vec& source,
437 const CellRangeType& cells,
440 auto s = source.begin();
441 auto c = cells.begin();
442 const auto e = cells.end();
443 for (; c != e; ++c, ++s) {
444 destination[*c] = *s;
455 #endif // OPM_INITSTATEEQUIL_HEADER_INCLUDED
void initStateEquil(const Grid &grid, const BlackoilPropertiesInterface &props, const Opm::DeckConstPtr deck, const Opm::EclipseStateConstPtr eclipseState, const double gravity, BlackoilState &state)
Definition: AnisotropicEikonal.hpp:43
Definition: BlackoilPhases.hpp:32
int phase_used[MaxNumPhases]
Definition: BlackoilPhases.hpp:39
virtual int numPhases() const =0
int numCells(const UnstructuredGrid &grid)
Get the number of cells of a grid.
Definition: BlackoilPropertiesInterface.hpp:37
int phase_pos[MaxNumPhases]
Definition: BlackoilPhases.hpp:40
const int * globalCell(const UnstructuredGrid &grid)
Get the local to global index mapping.
std::vector< std::vector< double > > phasePressures(const Grid &G, const Region ®, const CellRange &cells, const double grav)
Definition: initStateEquil_impl.hpp:594
virtual PhaseUsage phaseUsage() const =0
const double gravity
Definition: Units.hpp:120
std::vector< double > computeRs(const Grid &grid, const CellRangeType &cells, const std::vector< double > oil_pressure, const std::vector< double > &temperature, const Miscibility::RsFunction &rs_func, const std::vector< double > gas_saturation)
Definition: initStateEquil_impl.hpp:821
Definition: BlackoilPhases.hpp:32
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
std::vector< double > temperature(const Grid &, const Region &, const CellRange &cells)
Definition: initStateEquil_impl.hpp:667
std::vector< std::vector< double > > phaseSaturations(const Grid &G, const Region ®, const CellRange &cells, BlackoilPropertiesInterface &props, const std::vector< double > swat_init, std::vector< std::vector< double > > &phase_pressures)
Definition: initStateEquil_impl.hpp:677