23#ifndef OPM_INIT_STATE_EQUIL_IMPL_HPP
24#define OPM_INIT_STATE_EQUIL_IMPL_HPP
26#include <dune/grid/common/mcmgmapper.hh>
28#include <opm/common/OpmLog/OpmLog.hpp>
30#include <opm/grid/utility/RegionMapping.hpp>
32#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
33#include <opm/input/eclipse/EclipseState/Tables/PbvdTable.hpp>
34#include <opm/input/eclipse/EclipseState/Tables/PdvdTable.hpp>
35#include <opm/input/eclipse/EclipseState/Tables/RsconstTable.hpp>
36#include <opm/input/eclipse/EclipseState/Tables/RsvdTable.hpp>
37#include <opm/input/eclipse/EclipseState/Tables/RtempvdTable.hpp>
38#include <opm/input/eclipse/EclipseState/Tables/RvvdTable.hpp>
39#include <opm/input/eclipse/EclipseState/Tables/RvwvdTable.hpp>
40#include <opm/input/eclipse/EclipseState/Tables/SaltvdTable.hpp>
41#include <opm/input/eclipse/EclipseState/Tables/SaltpvdTable.hpp>
43#include <opm/input/eclipse/Units/UnitSystem.hpp>
45#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
46#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
53#include <fmt/format.h>
68template <
typename CellRange,
class Scalar>
70 const std::vector<std::pair<Scalar, Scalar>>&
cellZMinMax,
72 std::array<Scalar,2>& span)
74 span[0] = std::numeric_limits<Scalar>::max();
75 span[1] = std::numeric_limits<Scalar>::lowest();
85 for (
const auto& cell : cells) {
89 span[0] = comm.min(span[0]);
90 span[1] = comm.max(span[1]);
96 const int numIntervals,
97 std::vector<std::pair<Scalar, Scalar>>& subdiv)
99 const auto h = (right - left) / numIntervals;
102 for (
auto i = 0*numIntervals; i < numIntervals; ++i) {
103 const auto start = end;
104 end = left + (i + 1)*h;
106 subdiv.emplace_back((start + end) / 2, h);
110template <
typename CellID,
typename Scalar>
111std::vector<std::pair<Scalar, Scalar>>
113 const std::pair<Scalar, Scalar> topbot,
114 const int numIntervals)
116 auto subdiv = std::vector<std::pair<Scalar, Scalar>>{};
117 subdiv.reserve(2 * numIntervals);
119 if (topbot.first > topbot.second) {
120 throw std::out_of_range {
121 "Negative thickness (inverted top/bottom faces) in cell "
127 2*numIntervals, subdiv);
132template <
class Scalar,
class Element>
135 typedef typename Element::Geometry Geometry;
136 static constexpr int zCoord = Element::dimension - 1;
139 const Geometry& geometry = element.geometry();
140 const int corners = geometry.corners();
141 for (
int i=0; i < corners; ++i)
142 zz += geometry.corner(i)[zCoord];
147template <
class Scalar,
class Element>
150 typedef typename Element::Geometry Geometry;
151 static constexpr int xCoord = Element::dimension - 3;
152 static constexpr int yCoord = Element::dimension - 2;
157 const Geometry& geometry = element.geometry();
158 const int corners = geometry.corners();
159 for (
int i=0; i < corners; ++i) {
160 xx += geometry.corner(i)[xCoord];
161 yy += geometry.corner(i)[yCoord];
163 return std::make_pair(xx/corners, yy/corners);
166template <
class Scalar,
class Element>
167std::pair<Scalar,Scalar>
cellZSpan(
const Element& element)
169 typedef typename Element::Geometry Geometry;
170 static constexpr int zCoord = Element::dimension - 1;
174 const Geometry& geometry = element.geometry();
175 const int corners = geometry.corners();
176 assert(corners == 8);
177 for (
int i=0; i < 4; ++i)
178 bot += geometry.corner(i)[zCoord];
179 for (
int i=4; i < corners; ++i)
180 top += geometry.corner(i)[zCoord];
182 return std::make_pair(bot/4, top/4);
185template <
class Scalar,
class Element>
188 typedef typename Element::Geometry Geometry;
189 static constexpr int zCoord = Element::dimension - 1;
190 const Geometry& geometry = element.geometry();
191 const int corners = geometry.corners();
192 assert(corners == 8);
193 auto min = std::numeric_limits<Scalar>::max();
194 auto max = std::numeric_limits<Scalar>::lowest();
197 for (
int i=0; i < corners; ++i) {
198 min = std::min(min,
static_cast<Scalar
>(geometry.corner(i)[zCoord]));
199 max = std::max(max,
static_cast<Scalar
>(geometry.corner(i)[zCoord]));
201 return std::make_pair(min, max);
204template<
class Scalar>
206 Scalar& dipAngle, Scalar& dipAzimuth)
208 const auto& Xc = cellCorners.
X;
209 const auto& Yc = cellCorners.
Y;
210 const auto& Zc = cellCorners.
Z;
212 Scalar v1x = Xc[1] - Xc[0];
213 Scalar v1y = Yc[1] - Yc[0];
214 Scalar v1z = Zc[1] - Zc[0];
216 Scalar v2x = Xc[2] - Xc[0];
217 Scalar v2y = Yc[2] - Yc[0];
218 Scalar v2z = Zc[2] - Zc[0];
221 Scalar nx = v1y * v2z - v1z * v2y;
222 Scalar ny = v1z * v2x - v1x * v2z;
223 Scalar nz = v1x * v2y - v1y * v2x;
226 Scalar norm = std::hypot(nx, ny, nz);
234 dipAngle = std::acos(std::abs(nz));
237 if (std::abs(nx) > 1e-10 || std::abs(ny) > 1e-10) {
238 dipAzimuth = std::atan2(ny, nx);
240 dipAzimuth = std::fmod(dipAzimuth + 2*std::numbers::pi_v<Scalar>, 2*std::numbers::pi_v<Scalar>);
246 const Scalar maxDip = std::numbers::pi_v<Scalar>/2 -
static_cast<Scalar
>(1e-6);
247 dipAngle = std::min(dipAngle, maxDip);
255template <
class Scalar,
class Element>
258 typedef typename Element::Geometry Geometry;
259 const Geometry& geometry = element.geometry();
260 static constexpr int zCoord = Element::dimension - 1;
261 static constexpr int yCoord = Element::dimension - 2;
262 static constexpr int xCoord = Element::dimension - 3;
263 const int corners = geometry.corners();
264 assert(corners == 8);
265 std::array<Scalar, 8> X {};
266 std::array<Scalar, 8> Y {};
267 std::array<Scalar, 8> Z {};
269 for (
int i = 0; i < corners; ++i) {
270 auto corner = geometry.corner(i);
271 X[i] = corner[xCoord];
272 Y[i] = corner[yCoord];
273 Z[i] = corner[zCoord];
279template<
class Scalar>
281 Scalar dipAngle, Scalar dipAzimuth,
282 const std::array<Scalar, 3>& referencePoint)
289 Scalar dx = x - referencePoint[0];
290 Scalar dy = y - referencePoint[1];
291 Scalar dz = z - referencePoint[2];
294 if (std::abs(dipAngle) < 1e-10) {
295 return referencePoint[2] + dz;
299 Scalar pointAzimuth = std::atan2(dy, dx);
302 Scalar azimuthDiff = pointAzimuth - dipAzimuth;
305 Scalar lateralDist = std::hypot(dx, dy);
308 Scalar lateralInDipDir = lateralDist * std::cos(azimuthDiff);
313 Scalar tvd = referencePoint[2] + dz * std::cos(dipAngle) + lateralInDipDir * std::sin(dipAngle);
318template<
class Scalar,
class RHS>
320 const std::array<Scalar,2>& span,
326 const Scalar h = stepsize();
327 const Scalar h2 = h / 2;
328 const Scalar h6 = h / 6;
334 f_.push_back(f(span_[0], y0));
336 for (
int i = 0; i < N; ++i) {
337 const Scalar x = span_[0] + i*h;
338 const Scalar y = y_.back();
340 const Scalar k1 = f_[i];
341 const Scalar k2 = f(x + h2, y + h2*k1);
342 const Scalar k3 = f(x + h2, y + h2*k2);
343 const Scalar k4 = f(x + h, y + h*k3);
345 y_.push_back(y + h6*(k1 + 2*(k2 + k3) + k4));
346 f_.push_back(f(x + h, y_.back()));
349 assert (y_.size() ==
typename std::vector<Scalar>::size_type(N + 1));
352template<
class Scalar,
class RHS>
358 const Scalar h = stepsize();
359 int i = (x - span_[0]) / h;
360 const Scalar t = (x - (span_[0] + i*h)) / h;
363 if (i < 0) { i = 0; }
364 if (N_ <= i) { i = N_ - 1; }
366 const Scalar y0 = y_[i], y1 = y_[i + 1];
367 const Scalar f0 = f_[i], f1 = f_[i + 1];
369 Scalar u = (1 - 2*t) * (y1 - y0);
370 u += h * ((t - 1)*f0 + t*f1);
372 u += (1 - t)*y0 + t*y1;
377template<
class Scalar,
class RHS>
381 return (span_[1] - span_[0]) / N_;
384namespace PhasePressODE {
386template<
class Flu
idSystem>
388Water(
const TabulatedFunction& tempVdTable,
389 const TabulatedFunction& saltVdTable,
390 const int pvtRegionIdx,
391 const Scalar normGrav)
392 : tempVdTable_(tempVdTable)
393 , saltVdTable_(saltVdTable)
394 , pvtRegionIdx_(pvtRegionIdx)
399template<
class Flu
idSystem>
400typename Water<FluidSystem>::Scalar
403 const Scalar press)
const
405 return this->density(depth, press) * g_;
408template<
class Flu
idSystem>
409typename Water<FluidSystem>::Scalar
412 const Scalar press)
const
415 Scalar saltConcentration = saltVdTable_.eval(depth,
true);
416 Scalar temp = tempVdTable_.eval(depth,
true);
417 Scalar rho = FluidSystem::waterPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
422 rho *= FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
426template<
class Flu
idSystem,
class RS>
428Oil(
const TabulatedFunction& tempVdTable,
430 const int pvtRegionIdx,
431 const Scalar normGrav)
432 : tempVdTable_(tempVdTable)
434 , pvtRegionIdx_(pvtRegionIdx)
439template<
class Flu
idSystem,
class RS>
440typename Oil<FluidSystem,RS>::Scalar
443 const Scalar press)
const
445 return this->density(depth, press) * g_;
448template<
class Flu
idSystem,
class RS>
449typename Oil<FluidSystem,RS>::Scalar
452 const Scalar press)
const
454 const Scalar temp = tempVdTable_.eval(depth,
true);
456 if (FluidSystem::enableDissolvedGas() || FluidSystem::enableConstantRs())
457 rs = rs_(depth, press, temp);
460 if (rs >= FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press)) {
461 bOil = FluidSystem::oilPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
464 bOil = FluidSystem::oilPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rs);
466 Scalar rho = bOil * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
467 if (FluidSystem::enableDissolvedGas() || FluidSystem::enableConstantRs()) {
468 rho += rs * bOil * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
474template<
class Flu
idSystem,
class RV,
class RVW>
476Gas(
const TabulatedFunction& tempVdTable,
479 const int pvtRegionIdx,
480 const Scalar normGrav)
481 : tempVdTable_(tempVdTable)
484 , pvtRegionIdx_(pvtRegionIdx)
489template<
class Flu
idSystem,
class RV,
class RVW>
490typename Gas<FluidSystem,RV,RVW>::Scalar
493 const Scalar press)
const
495 return this->density(depth, press) * g_;
498template<
class Flu
idSystem,
class RV,
class RVW>
499typename Gas<FluidSystem,RV,RVW>::Scalar
502 const Scalar press)
const
504 const Scalar temp = tempVdTable_.eval(depth,
true);
506 if (FluidSystem::enableVaporizedOil())
507 rv = rv_(depth, press, temp);
510 if (FluidSystem::enableVaporizedWater())
511 rvw = rvw_(depth, press, temp);
515 if (FluidSystem::enableVaporizedOil() && FluidSystem::enableVaporizedWater()) {
516 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)
517 && rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press))
519 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
521 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rv, rvw);
523 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
524 rho += rv * bGas * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_)
525 + rvw * bGas * FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
529 if (FluidSystem::enableVaporizedOil()){
530 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)) {
531 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
533 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
539 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
540 rho += rv * bGas * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
544 if (FluidSystem::enableVaporizedWater()){
545 if (rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press)) {
546 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
549 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
555 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
556 rho += rvw * bGas * FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
561 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp,
565 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
572template<
class Flu
idSystem,
class Region>
574PressureTable<FluidSystem,Region>::
575PressureFunction<ODE>::PressureFunction(
const ODE& ode,
581 this->value_[Direction::Up] = std::make_unique<Distribution>
582 (ode, VSpan {{ ic.depth, span[0] }}, ic.pressure, nsample);
584 this->value_[Direction::Down] = std::make_unique<Distribution>
585 (ode, VSpan {{ ic.depth, span[1] }}, ic.pressure, nsample);
588template<
class Flu
idSystem,
class Region>
590PressureTable<FluidSystem,Region>::
591PressureFunction<ODE>::PressureFunction(
const PressureFunction& rhs)
592 : initial_(rhs.initial_)
594 this->value_[Direction::Up] =
595 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
597 this->value_[Direction::Down] =
598 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
601template<
class Flu
idSystem,
class Region>
603typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
608 this->initial_ = rhs.initial_;
610 this->value_[Direction::Up] =
611 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
613 this->value_[Direction::Down] =
614 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
619template<
class Flu
idSystem,
class Region>
621typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
626 this->initial_ = rhs.initial_;
627 this->value_ = std::move(rhs.value_);
632template<
class Flu
idSystem,
class Region>
635PressureTable<FluidSystem,Region>::
636PressureFunction<ODE>::
637value(
const Scalar depth)
const
639 if (depth < this->initial_.depth) {
641 return (*this->value_[Direction::Up])(depth);
643 else if (depth > this->initial_.depth) {
645 return (*this->value_[Direction::Down])(depth);
649 return this->initial_.pressure;
654template<
class Flu
idSystem,
class Region>
655template<
typename PressFunc>
656void PressureTable<FluidSystem,Region>::
657checkPtr(
const PressFunc* phasePress,
658 const std::string& phaseName)
const
660 if (phasePress !=
nullptr) {
return; }
662 throw std::invalid_argument {
663 "Phase pressure function for \"" + phaseName
664 +
"\" most not be null"
668template<
class Flu
idSystem,
class Region>
669typename PressureTable<FluidSystem,Region>::Strategy
670PressureTable<FluidSystem,Region>::
671selectEquilibrationStrategy(
const Region& reg)
const
673 if (!this->oilActive()) {
674 if (reg.datum() > reg.zwoc()) {
675 return &PressureTable::equil_WOG;
677 return &PressureTable::equil_GOW;
680 if (reg.datum() > reg.zwoc()) {
681 return &PressureTable::equil_WOG;
683 else if (reg.datum() < reg.zgoc()) {
684 return &PressureTable::equil_GOW;
687 return &PressureTable::equil_OWG;
691template<
class Flu
idSystem,
class Region>
692void PressureTable<FluidSystem,Region>::
693copyInPointers(
const PressureTable& rhs)
695 if (rhs.oil_ !=
nullptr) {
696 this->oil_ = std::make_unique<OPress>(*rhs.oil_);
699 if (rhs.gas_ !=
nullptr) {
700 this->gas_ = std::make_unique<GPress>(*rhs.gas_);
703 if (rhs.wat_ !=
nullptr) {
704 this->wat_ = std::make_unique<WPress>(*rhs.wat_);
708template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
711 const std::vector<Scalar>& swatInit)
712 : matLawMgr_(matLawMgr)
713 , swatInit_ (swatInit)
717template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
720 : matLawMgr_(rhs.matLawMgr_)
721 , swatInit_ (rhs.swatInit_)
723 , press_ (rhs.press_)
726 this->setEvaluationPoint(*rhs.evalPt_.position,
728 *rhs.evalPt_.ptable);
731template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
738 this->setEvaluationPoint(x, reg, ptable);
739 this->initializePhaseQuantities();
741 if (ptable.
gasActive()) { this->deriveGasSat(); }
743 if (ptable.
waterActive()) { this->deriveWaterSat(); }
746 if (this->isOverlappingTransition()) {
747 this->fixUnphysicalTransition();
750 if (ptable.
oilActive()) { this->deriveOilSat(); }
752 this->accountForScaledSaturations();
757template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
761 const PTable& ptable)
763 this->evalPt_.position = &x;
764 this->evalPt_.region = ®
765 this->evalPt_.ptable = &ptable;
768template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
769void PhaseSaturations<MaterialLawManager,FluidSystem,Region,CellID>::
770initializePhaseQuantities()
773 this->press_.reset();
775 const auto depth = this->evalPt_.position->depth;
776 const auto& ptable = *this->evalPt_.ptable;
778 if (ptable.oilActive()) {
779 this->press_.oil = ptable.oil(depth);
782 if (ptable.gasActive()) {
783 this->press_.gas = ptable.gas(depth);
786 if (ptable.waterActive()) {
787 this->press_.water = ptable.water(depth);
791template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
792void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveOilSat()
794 this->sat_.oil = 1.0 - this->sat_.water - this->sat_.gas;
797template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
798void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveGasSat()
800 auto& sg = this->sat_.gas;
802 const auto isIncr =
true;
803 const auto oilActive = this->evalPt_.ptable->oilActive();
805 if (this->isConstCapPress(this->gasPos())) {
809 const auto gas_contact = oilActive? this->evalPt_.region->zgoc() : this->evalPt_.region->zwoc();
810 sg = this->fromDepthTable(gas_contact,
811 this->gasPos(), isIncr);
821 const auto pw = oilActive? this->press_.oil : this->press_.water;
822 const auto pcgo = this->press_.gas - pw;
823 sg = this->invertCapPress(pcgo, this->gasPos(), isIncr);
827template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
828void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveWaterSat()
830 auto& sw = this->sat_.water;
832 const auto oilActive = this->evalPt_.ptable->oilActive();
835 sw = 1.0 - this->sat_.gas;
838 const auto isIncr =
false;
840 if (this->isConstCapPress(this->waterPos())) {
844 sw = this->fromDepthTable(this->evalPt_.region->zwoc(),
845 this->waterPos(), isIncr);
857 const auto pcow = this->press_.oil - this->press_.water;
859 if (this->swatInit_.empty()) {
860 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
863 auto [swout, newSwatInit] = this->applySwatInit(pcow);
865 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
874template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
875void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
876fixUnphysicalTransition()
878 auto& sg = this->sat_.gas;
879 auto& sw = this->sat_.water;
887 const auto pcgw = this->press_.gas - this->press_.water;
888 if (! this->swatInit_.empty()) {
892 auto [swout, newSwatInit] = this->applySwatInit(pcgw, sw);
894 const auto isIncr =
false;
895 sw = this->invertCapPress(pcgw, this->waterPos(), isIncr);
902 sw = satFromSumOfPcs<FluidSystem>
903 (this->matLawMgr_, this->waterPos(), this->gasPos(),
904 this->evalPt_.position->cell, pcgw);
907 this->fluidState_.setSaturation(this->oilPos(), 1.0 - sw - sg);
908 this->fluidState_.setSaturation(this->gasPos(), sg);
909 this->fluidState_.setSaturation(this->waterPos(), this->evalPt_
910 .ptable->waterActive() ? sw : 0.0);
913 this->computeMaterialLawCapPress();
914 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
917template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
918void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
919accountForScaledSaturations()
921 const auto gasActive = this->evalPt_.ptable->gasActive();
922 const auto watActive = this->evalPt_.ptable->waterActive();
923 const auto oilActive = this->evalPt_.ptable->oilActive();
925 auto sg = gasActive? this->sat_.gas : 0.0;
926 auto sw = watActive? this->sat_.water : 0.0;
927 auto so = oilActive? this->sat_.oil : 0.0;
929 this->fluidState_.setSaturation(this->waterPos(), sw);
930 this->fluidState_.setSaturation(this->oilPos(), so);
931 this->fluidState_.setSaturation(this->gasPos(), sg);
933 const auto& scaledDrainageInfo = this->matLawMgr_
934 .oilWaterScaledEpsInfoDrainage(this->evalPt_.position->cell);
936 const auto thresholdSat = 1.0e-6;
937 if (watActive && ((sw + thresholdSat) > scaledDrainageInfo.Swu)) {
941 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swu);
943 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swu);
944 }
else if (gasActive) {
945 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swu);
947 sw = scaledDrainageInfo.Swu;
948 this->computeMaterialLawCapPress();
952 this->press_.oil = this->press_.water + this->materialLawCapPressOilWater();
955 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
959 if (gasActive && ((sg + thresholdSat) > scaledDrainageInfo.Sgu)) {
963 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgu);
965 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgu);
966 }
else if (watActive) {
967 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgu);
969 sg = scaledDrainageInfo.Sgu;
970 this->computeMaterialLawCapPress();
974 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
977 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
981 if (watActive && ((sw - thresholdSat) < scaledDrainageInfo.Swl)) {
985 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swl);
987 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swl);
988 }
else if (gasActive) {
989 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swl);
991 sw = scaledDrainageInfo.Swl;
992 this->computeMaterialLawCapPress();
996 this->press_.water = this->press_.oil - this->materialLawCapPressOilWater();
999 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
1003 if (gasActive && ((sg - thresholdSat) < scaledDrainageInfo.Sgl)) {
1007 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgl);
1009 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgl);
1010 }
else if (watActive) {
1011 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgl);
1013 sg = scaledDrainageInfo.Sgl;
1014 this->computeMaterialLawCapPress();
1018 this->press_.gas = this->press_.oil + this->materialLawCapPressGasOil();
1021 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
1026template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1027std::pair<typename FluidSystem::Scalar, bool>
1028PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1029applySwatInit(
const Scalar pcow)
1031 return this->applySwatInit(pcow, this->swatInit_[this->evalPt_.position->cell]);
1034template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1035std::pair<typename FluidSystem::Scalar, bool>
1036PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1037applySwatInit(
const Scalar pcow,
const Scalar sw)
1039 return this->matLawMgr_.applySwatinit(this->evalPt_.position->cell, pcow, sw);
1042template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1043void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1044computeMaterialLawCapPress()
1046 const auto& matParams = this->matLawMgr_
1047 .materialLawParams(this->evalPt_.position->cell);
1049 this->matLawCapPress_.fill(0.0);
1050 MaterialLaw::capillaryPressures(this->matLawCapPress_,
1051 matParams, this->fluidState_);
1054template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1055typename FluidSystem::Scalar
1056PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1057materialLawCapPressGasOil()
const
1059 return this->matLawCapPress_[this->oilPos()]
1060 + this->matLawCapPress_[this->gasPos()];
1063template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1064typename FluidSystem::Scalar
1065PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1066materialLawCapPressOilWater()
const
1068 return this->matLawCapPress_[this->oilPos()]
1069 - this->matLawCapPress_[this->waterPos()];
1072template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1073typename FluidSystem::Scalar
1074PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1075materialLawCapPressGasWater()
const
1077 return this->matLawCapPress_[this->gasPos()]
1078 - this->matLawCapPress_[this->waterPos()];
1081template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1082bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1083isConstCapPress(
const PhaseIdx phaseIdx)
const
1085 return isConstPc<FluidSystem>
1086 (this->matLawMgr_, phaseIdx, this->evalPt_.position->cell);
1089template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1090bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1091isOverlappingTransition()
const
1093 return this->evalPt_.ptable->gasActive()
1094 && this->evalPt_.ptable->waterActive()
1095 && ((this->sat_.gas + this->sat_.water) > 1.0);
1098template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1099typename FluidSystem::Scalar
1100PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1101fromDepthTable(
const Scalar contactdepth,
1102 const PhaseIdx phasePos,
1103 const bool isincr)
const
1105 return satFromDepth<FluidSystem>
1106 (this->matLawMgr_, this->evalPt_.position->depth,
1107 contactdepth,
static_cast<int>(phasePos),
1108 this->evalPt_.position->cell, isincr);
1111template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1112typename FluidSystem::Scalar
1113PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1114invertCapPress(
const Scalar pc,
1115 const PhaseIdx phasePos,
1116 const bool isincr)
const
1118 return satFromPc<FluidSystem>
1119 (this->matLawMgr_,
static_cast<int>(phasePos),
1120 this->evalPt_.position->cell, pc, isincr);
1123template<
class Flu
idSystem,
class Region>
1126 const int samplePoints)
1128 , nsample_(samplePoints)
1132template <
class Flu
idSystem,
class Region>
1135 : gravity_(rhs.gravity_)
1136 , nsample_(rhs.nsample_)
1138 this->copyInPointers(rhs);
1141template <
class Flu
idSystem,
class Region>
1144 : gravity_(rhs.gravity_)
1145 , nsample_(rhs.nsample_)
1146 , oil_ (std::move(rhs.oil_))
1147 , gas_ (std::move(rhs.gas_))
1148 , wat_ (std::move(rhs.wat_))
1152template <
class Flu
idSystem,
class Region>
1157 this->gravity_ = rhs.gravity_;
1158 this->nsample_ = rhs.nsample_;
1159 this->copyInPointers(rhs);
1164template <
class Flu
idSystem,
class Region>
1169 this->gravity_ = rhs.gravity_;
1170 this->nsample_ = rhs.nsample_;
1172 this->oil_ = std::move(rhs.oil_);
1173 this->gas_ = std::move(rhs.gas_);
1174 this->wat_ = std::move(rhs.wat_);
1179template <
class Flu
idSystem,
class Region>
1185 auto equil = this->selectEquilibrationStrategy(reg);
1187 (this->*equil)(reg, span);
1190template <
class Flu
idSystem,
class Region>
1194 return FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1197template <
class Flu
idSystem,
class Region>
1201 return FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
1204template <
class Flu
idSystem,
class Region>
1208 return FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1211template <
class Flu
idSystem,
class Region>
1212typename FluidSystem::Scalar
1216 this->checkPtr(this->oil_.get(),
"OIL");
1218 return this->oil_->value(depth);
1221template <
class Flu
idSystem,
class Region>
1222typename FluidSystem::Scalar
1226 this->checkPtr(this->gas_.get(),
"GAS");
1228 return this->gas_->value(depth);
1232template <
class Flu
idSystem,
class Region>
1233typename FluidSystem::Scalar
1237 this->checkPtr(this->wat_.get(),
"WATER");
1239 return this->wat_->value(depth);
1242template <
class Flu
idSystem,
class Region>
1244equil_WOG(
const Region& reg,
const VSpan& span)
1249 if (! this->waterActive()) {
1250 throw std::invalid_argument {
1251 "Don't know how to interpret EQUIL datum depth in "
1252 "WATER zone in model without active water phase"
1257 const auto ic =
typename WPress::InitCond {
1258 reg.datum(), reg.pressure()
1261 this->makeWatPressure(ic, reg, span);
1264 if (this->oilActive()) {
1266 const auto ic =
typename OPress::InitCond {
1268 this->water(reg.zwoc()) + reg.pcowWoc()
1271 this->makeOilPressure(ic, reg, span);
1274 if (this->gasActive() && this->oilActive()) {
1276 const auto ic =
typename GPress::InitCond {
1278 this->oil(reg.zgoc()) + reg.pcgoGoc()
1281 this->makeGasPressure(ic, reg, span);
1282 }
else if (this->gasActive() && !this->oilActive()) {
1284 const auto ic =
typename GPress::InitCond {
1286 this->water(reg.zwoc()) + reg.pcowWoc()
1288 this->makeGasPressure(ic, reg, span);
1292template <
class Flu
idSystem,
class Region>
1293void PressureTable<FluidSystem, Region>::
1294equil_GOW(
const Region& reg,
const VSpan& span)
1299 if (! this->gasActive()) {
1300 throw std::invalid_argument {
1301 "Don't know how to interpret EQUIL datum depth in "
1302 "GAS zone in model without active gas phase"
1307 const auto ic =
typename GPress::InitCond {
1308 reg.datum(), reg.pressure()
1311 this->makeGasPressure(ic, reg, span);
1314 if (this->oilActive()) {
1316 const auto ic =
typename OPress::InitCond {
1318 this->gas(reg.zgoc()) - reg.pcgoGoc()
1320 this->makeOilPressure(ic, reg, span);
1323 if (this->waterActive() && this->oilActive()) {
1325 const auto ic =
typename WPress::InitCond {
1327 this->oil(reg.zwoc()) - reg.pcowWoc()
1330 this->makeWatPressure(ic, reg, span);
1331 }
else if (this->waterActive() && !this->oilActive()) {
1333 const auto ic =
typename WPress::InitCond {
1335 this->gas(reg.zwoc()) - reg.pcowWoc()
1337 this->makeWatPressure(ic, reg, span);
1341template <
class Flu
idSystem,
class Region>
1342void PressureTable<FluidSystem, Region>::
1343equil_OWG(
const Region& reg,
const VSpan& span)
1348 if (! this->oilActive()) {
1349 throw std::invalid_argument {
1350 "Don't know how to interpret EQUIL datum depth in "
1351 "OIL zone in model without active oil phase"
1356 const auto ic =
typename OPress::InitCond {
1357 reg.datum(), reg.pressure()
1360 this->makeOilPressure(ic, reg, span);
1363 if (this->waterActive()) {
1365 const auto ic =
typename WPress::InitCond {
1367 this->oil(reg.zwoc()) - reg.pcowWoc()
1370 this->makeWatPressure(ic, reg, span);
1373 if (this->gasActive()) {
1375 const auto ic =
typename GPress::InitCond {
1377 this->oil(reg.zgoc()) + reg.pcgoGoc()
1379 this->makeGasPressure(ic, reg, span);
1383template <
class Flu
idSystem,
class Region>
1384void PressureTable<FluidSystem, Region>::
1385makeOilPressure(
const typename OPress::InitCond& ic,
1389 const auto drho = OilPressODE {
1390 reg.tempVdTable(), reg.dissolutionCalculator(),
1391 reg.pvtIdx(), this->gravity_
1394 this->oil_ = std::make_unique<OPress>(drho, ic, this->nsample_, span);
1397template <
class Flu
idSystem,
class Region>
1398void PressureTable<FluidSystem, Region>::
1399makeGasPressure(
const typename GPress::InitCond& ic,
1403 const auto drho = GasPressODE {
1404 reg.tempVdTable(), reg.evaporationCalculator(), reg.waterEvaporationCalculator(),
1405 reg.pvtIdx(), this->gravity_
1408 this->gas_ = std::make_unique<GPress>(drho, ic, this->nsample_, span);
1411template <
class Flu
idSystem,
class Region>
1412void PressureTable<FluidSystem, Region>::
1413makeWatPressure(
const typename WPress::InitCond& ic,
1417 const auto drho = WatPressODE {
1418 reg.tempVdTable(), reg.saltVdTable(), reg.pvtIdx(), this->gravity_
1421 this->wat_ = std::make_unique<WPress>(drho, ic, this->nsample_, span);
1426namespace DeckDependent {
1428std::vector<EquilRecord>
1431 const auto& init = state.getInitConfig();
1433 if(!init.hasEquil()) {
1434 throw std::domain_error(
"Deck does not provide equilibration data.");
1437 const auto& equil = init.getEquil();
1438 return { equil.begin(), equil.end() };
1441template<
class Gr
idView>
1444 const GridView& gridview)
1446 std::vector<int> eqlnum(gridview.size(0), 0);
1448 if (eclipseState.fieldProps().has_int(
"EQLNUM")) {
1449 const auto& e = eclipseState.fieldProps().get_int(
"EQLNUM");
1450 std::ranges::transform(e, eqlnum.begin(), [](
int n) { return n - 1; });
1453 const int num_regions = eclipseState.getTableManager().getEqldims().getNumEquilRegions();
1454 if (std::ranges::any_of(eqlnum, [num_regions](
int n){
return n >= num_regions;})) {
1455 throw std::runtime_error(
"Values larger than maximum Equil regions " +
1458 if (std::ranges::any_of(eqlnum, [](
int n){
return n < 0;})) {
1459 throw std::runtime_error(
"zero or negative values provided in EQLNUM");
1466template<
class FluidSystem,
1469 class ElementMapper,
1470 class CartesianIndexMapper>
1471template<
class MaterialLawManager>
1472InitialStateComputer<FluidSystem,
1476 CartesianIndexMapper>::
1477InitialStateComputer(MaterialLawManager& materialLawManager,
1478 const EclipseState& eclipseState,
1480 const GridView& gridView,
1481 const CartesianIndexMapper& cartMapper,
1483 const int num_pressure_points,
1484 const bool applySwatInit)
1485 : temperature_(grid.size(0), eclipseState.getTableManager().rtemp()),
1486 saltConcentration_(grid.size(0)),
1487 saltSaturation_(grid.size(0)),
1488 pp_(FluidSystem::numPhases,
1489 std::vector<Scalar>(grid.size(0))),
1490 sat_(FluidSystem::numPhases,
1491 std::vector<Scalar>(grid.size(0))),
1495 cartesianIndexMapper_(cartMapper),
1496 num_pressure_points_(num_pressure_points)
1499 if (applySwatInit) {
1500 if (eclipseState.fieldProps().has_double(
"SWATINIT")) {
1501 if constexpr (std::is_same_v<Scalar,double>) {
1502 swatInit_ = eclipseState.fieldProps().get_double(
"SWATINIT");
1504 const auto& input = eclipseState.fieldProps().get_double(
"SWATINIT");
1505 swatInit_.resize(input.size());
1506 std::ranges::copy(input, swatInit_.begin());
1513 const auto& num_aquifers = eclipseState.aquifer().numericalAquifers();
1514 updateCellProps_(gridView, num_aquifers);
1517 const std::vector<EquilRecord> rec =
getEquil(eclipseState);
1518 const auto& tables = eclipseState.getTableManager();
1520 const RegionMapping<> eqlmap(
equilnum(eclipseState, grid));
1521 const int invalidRegion = -1;
1522 regionPvtIdx_.resize(rec.size(), invalidRegion);
1523 setRegionPvtIdx(eclipseState, eqlmap);
1526 rsFunc_.reserve(rec.size());
1528 auto getArray = [](
const std::vector<double>& input)
1530 if constexpr (std::is_same_v<Scalar,double>) {
1533 std::vector<Scalar> output;
1534 output.resize(input.size());
1535 std::ranges::copy(input, output.begin());
1540 if (FluidSystem::enableDissolvedGas()) {
1541 for (std::size_t i = 0; i < rec.size(); ++i) {
1542 if (eqlmap.cells(i).empty()) {
1546 const int pvtIdx = regionPvtIdx_[i];
1547 if (!rec[i].liveOilInitConstantRs()) {
1548 const TableContainer& rsvdTables = tables.getRsvdTables();
1549 const TableContainer& pbvdTables = tables.getPbvdTables();
1550 if (rsvdTables.size() > 0) {
1551 const RsvdTable& rsvdTable = rsvdTables.getTable<RsvdTable>(i);
1552 auto depthColumn = getArray(rsvdTable.getColumn(
"DEPTH").vectorCopy());
1553 auto rsColumn = getArray(rsvdTable.getColumn(
"RS").vectorCopy());
1555 depthColumn, rsColumn));
1556 }
else if (pbvdTables.size() > 0) {
1557 const PbvdTable& pbvdTable = pbvdTables.getTable<PbvdTable>(i);
1558 auto depthColumn = getArray(pbvdTable.getColumn(
"DEPTH").vectorCopy());
1559 auto pbubColumn = getArray(pbvdTable.getColumn(
"PBUB").vectorCopy());
1561 depthColumn, pbubColumn));
1564 throw std::runtime_error(
"Cannot initialise: RSVD or PBVD table not available.");
1569 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1570 throw std::runtime_error(
"Cannot initialise: when no explicit RSVD table is given, \n"
1571 "datum depth must be at the gas-oil-contact. "
1572 "In EQUIL region "+
std::to_string(i + 1)+
" (counting from 1), this does not hold.");
1574 const Scalar pContact = rec[i].datumDepthPressure();
1575 const Scalar TContact = 273.15 + 20;
1580 else if (FluidSystem::enableConstantRs() && tables.hasTables(
"RSCONST")) {
1581 const auto& rsconstTables = tables.getRsconstTables();
1583 if (rsconstTables.empty()) {
1584 for (std::size_t i = 0; i < rec.size(); ++i) {
1590 const auto& rsconstTable = rsconstTables.getTable<RsconstTable>(0);
1592 const auto rsConst = rsconstTable.getRsColumn().front();
1593 const auto pBub = rsconstTable.getPbubColumn().front();
1595 const auto& units = eclipseState.getUnits();
1597 OpmLog::info(fmt::format(
"Using RSCONST keyword: Rs = {:.2} [{}], Pb = {:.2} [{}]",
1598 units.from_si(UnitSystem::measure::gas_oil_ratio, rsConst),
1599 units.name (UnitSystem::measure::gas_oil_ratio),
1600 units.from_si(UnitSystem::measure::pressure, pBub),
1601 units.name (UnitSystem::measure::pressure)));
1603 for (std::size_t i = 0; i < rec.size(); ++i) {
1609 for (std::size_t i = 0; i < rec.size(); ++i) {
1615 rvFunc_.reserve(rec.size());
1616 if (FluidSystem::enableVaporizedOil()) {
1617 for (std::size_t i = 0; i < rec.size(); ++i) {
1618 if (eqlmap.cells(i).empty()) {
1622 const int pvtIdx = regionPvtIdx_[i];
1623 if (!rec[i].wetGasInitConstantRv()) {
1624 const TableContainer& rvvdTables = tables.getRvvdTables();
1625 const TableContainer& pdvdTables = tables.getPdvdTables();
1627 if (rvvdTables.size() > 0) {
1628 const RvvdTable& rvvdTable = rvvdTables.getTable<RvvdTable>(i);
1629 auto depthColumn = getArray(rvvdTable.getColumn(
"DEPTH").vectorCopy());
1630 auto rvColumn = getArray(rvvdTable.getColumn(
"RV").vectorCopy());
1632 depthColumn, rvColumn));
1633 }
else if (pdvdTables.size() > 0) {
1634 const PdvdTable& pdvdTable = pdvdTables.getTable<PdvdTable>(i);
1635 auto depthColumn = getArray(pdvdTable.getColumn(
"DEPTH").vectorCopy());
1636 auto pdewColumn = getArray(pdvdTable.getColumn(
"PDEW").vectorCopy());
1638 depthColumn, pdewColumn));
1640 throw std::runtime_error(
"Cannot initialise: RVVD or PDCD table not available.");
1644 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1645 throw std::runtime_error(
1646 "Cannot initialise: when no explicit RVVD table is given, \n"
1647 "datum depth must be at the gas-oil-contact. "
1648 "In EQUIL region "+
std::to_string(i + 1)+
" (counting from 1), this does not hold.");
1650 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1651 const Scalar TContact = 273.15 + 20;
1657 for (std::size_t i = 0; i < rec.size(); ++i) {
1662 rvwFunc_.reserve(rec.size());
1663 if (FluidSystem::enableVaporizedWater()) {
1664 for (std::size_t i = 0; i < rec.size(); ++i) {
1665 if (eqlmap.cells(i).empty()) {
1669 const int pvtIdx = regionPvtIdx_[i];
1670 if (!rec[i].humidGasInitConstantRvw()) {
1671 const TableContainer& rvwvdTables = tables.getRvwvdTables();
1673 if (rvwvdTables.size() > 0) {
1674 const RvwvdTable& rvwvdTable = rvwvdTables.getTable<RvwvdTable>(i);
1675 auto depthColumn = getArray(rvwvdTable.getColumn(
"DEPTH").vectorCopy());
1676 auto rvwvdColumn = getArray(rvwvdTable.getColumn(
"RVWVD").vectorCopy());
1678 depthColumn, rvwvdColumn));
1680 throw std::runtime_error(
"Cannot initialise: RVWVD table not available.");
1684 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1686 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1688 const auto msg =
"No explicit RVWVD table is given for EQUIL region " +
std::to_string(i + 1) +
". \n"
1689 "and datum depth is not at the gas-oil-contact. \n"
1690 "Rvw is set to 0.0 in all cells. \n";
1691 OpmLog::warning(msg);
1695 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1696 const Scalar TContact = 273.15 + 20;
1703 if (rec[i].waterOilContactDepth() != rec[i].datumDepth()) {
1705 const auto msg =
"No explicit RVWVD table is given for EQUIL region " +
std::to_string(i + 1) +
". \n"
1706 "and datum depth is not at the gas-water-contact. \n"
1707 "Rvw is set to 0.0 in all cells. \n";
1708 OpmLog::warning(msg);
1711 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].waterOilContactCapillaryPressure();
1712 const Scalar TContact = 273.15 + 20;
1720 for (std::size_t i = 0; i < rec.size(); ++i) {
1726 updateInitialTemperature_(eclipseState, eqlmap);
1729 updateInitialSaltConcentration_(eclipseState, eqlmap);
1732 updateInitialSaltSaturation_(eclipseState, eqlmap);
1735 const auto& comm = grid.comm();
1736 calcPressSatRsRv(eqlmap, rec, materialLawManager, gridView, comm, grav);
1739 applyNumericalAquifers_(gridView, num_aquifers,
1740 eclipseState.runspec().co2Storage() ||
1741 eclipseState.runspec().h2Storage());
1747template<
class FluidSystem,
1750 class ElementMapper,
1751 class CartesianIndexMapper>
1757 CartesianIndexMapper>::
1758updateInitialTemperature_(
const EclipseState& eclState,
const RMap& reg)
1760 const int numEquilReg = rsFunc_.size();
1761 tempVdTable_.resize(numEquilReg);
1762 const auto& tables = eclState.getTableManager();
1763 if (!tables.hasTables(
"RTEMPVD")) {
1764 std::vector<Scalar> x = {0.0,1.0};
1765 std::vector<Scalar> y = {
static_cast<Scalar
>(tables.rtemp()),
1766 static_cast<Scalar
>(tables.rtemp())};
1767 for (
auto& table : this->tempVdTable_) {
1768 table.setXYContainers(x, y);
1771 const TableContainer& tempvdTables = tables.getRtempvdTables();
1772 for (std::size_t i = 0; i < tempvdTables.size(); ++i) {
1773 const RtempvdTable& tempvdTable = tempvdTables.getTable<RtempvdTable>(i);
1774 tempVdTable_[i].setXYContainers(tempvdTable.getDepthColumn(), tempvdTable.getTemperatureColumn());
1775 const auto& cells = reg.cells(i);
1776 for (
const auto& cell : cells) {
1777 const Scalar depth = cellCenterDepth_[cell];
1778 this->temperature_[cell] = tempVdTable_[i].eval(depth,
true);
1784template<
class FluidSystem,
1787 class ElementMapper,
1788 class CartesianIndexMapper>
1790void InitialStateComputer<FluidSystem,
1794 CartesianIndexMapper>::
1795updateInitialSaltConcentration_(
const EclipseState& eclState,
const RMap& reg)
1797 const int numEquilReg = rsFunc_.size();
1798 saltVdTable_.resize(numEquilReg);
1799 const auto& tables = eclState.getTableManager();
1800 const TableContainer& saltvdTables = tables.getSaltvdTables();
1803 if (saltvdTables.empty()) {
1804 std::vector<Scalar> x = {0.0,1.0};
1805 std::vector<Scalar> y = {0.0,0.0};
1806 for (
auto& table : this->saltVdTable_) {
1807 table.setXYContainers(x, y);
1810 for (std::size_t i = 0; i < saltvdTables.size(); ++i) {
1811 const SaltvdTable& saltvdTable = saltvdTables.getTable<SaltvdTable>(i);
1812 saltVdTable_[i].setXYContainers(saltvdTable.getDepthColumn(), saltvdTable.getSaltColumn());
1814 const auto& cells = reg.cells(i);
1815 for (
const auto& cell : cells) {
1816 const Scalar depth = cellCenterDepth_[cell];
1817 this->saltConcentration_[cell] = saltVdTable_[i].eval(depth,
true);
1823template<
class FluidSystem,
1826 class ElementMapper,
1827 class CartesianIndexMapper>
1829void InitialStateComputer<FluidSystem,
1833 CartesianIndexMapper>::
1834updateInitialSaltSaturation_(
const EclipseState& eclState,
const RMap& reg)
1836 const int numEquilReg = rsFunc_.size();
1837 saltpVdTable_.resize(numEquilReg);
1838 const auto& tables = eclState.getTableManager();
1839 const TableContainer& saltpvdTables = tables.getSaltpvdTables();
1841 for (std::size_t i = 0; i < saltpvdTables.size(); ++i) {
1842 const SaltpvdTable& saltpvdTable = saltpvdTables.getTable<SaltpvdTable>(i);
1843 saltpVdTable_[i].setXYContainers(saltpvdTable.getDepthColumn(), saltpvdTable.getSaltpColumn());
1845 const auto& cells = reg.cells(i);
1846 for (
const auto& cell : cells) {
1847 const Scalar depth = cellCenterDepth_[cell];
1848 this->saltSaturation_[cell] = saltpVdTable_[i].eval(depth,
true);
1853template<
class FluidSystem,
1856 class ElementMapper,
1857 class CartesianIndexMapper>
1858void InitialStateComputer<FluidSystem,
1862 CartesianIndexMapper>::
1863updateCellProps_(
const GridView& gridView,
1864 const NumericalAquifers& aquifer)
1866 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1867 int numElements = gridView.size(0);
1868 cellCenterDepth_.resize(numElements);
1869 cellCenterXY_.resize(numElements);
1870 cellCorners_.resize(numElements);
1871 cellZSpan_.resize(numElements);
1872 cellZMinMax_.resize(numElements);
1874 auto elemIt = gridView.template begin<0>();
1875 const auto& elemEndIt = gridView.template end<0>();
1876 const auto num_aqu_cells = aquifer.allAquiferCells();
1877 for (; elemIt != elemEndIt; ++elemIt) {
1878 const Element& element = *elemIt;
1879 const unsigned int elemIdx = elemMapper.index(element);
1880 cellCenterDepth_[elemIdx] = Details::cellCenterDepth<Scalar>(element);
1881 cellCenterXY_[elemIdx] = Details::cellCenterXY<Scalar>(element);
1882 cellCorners_[elemIdx] = Details::getCellCornerXY<Scalar>(element);
1883 const auto cartIx = cartesianIndexMapper_.cartesianIndex(elemIdx);
1884 cellZSpan_[elemIdx] = Details::cellZSpan<Scalar>(element);
1885 cellZMinMax_[elemIdx] = Details::cellZMinMax<Scalar>(element);
1886 if (!num_aqu_cells.empty()) {
1887 const auto search = num_aqu_cells.find(cartIx);
1888 if (search != num_aqu_cells.end()) {
1889 const auto* aqu_cell = num_aqu_cells.at(cartIx);
1890 const Scalar depth_change_num_aqu = aqu_cell->depth - cellCenterDepth_[elemIdx];
1891 cellCenterDepth_[elemIdx] += depth_change_num_aqu;
1892 cellZSpan_[elemIdx].first += depth_change_num_aqu;
1893 cellZSpan_[elemIdx].second += depth_change_num_aqu;
1894 cellZMinMax_[elemIdx].first += depth_change_num_aqu;
1895 cellZMinMax_[elemIdx].second += depth_change_num_aqu;
1901template<
class FluidSystem,
1904 class ElementMapper,
1905 class CartesianIndexMapper>
1906void InitialStateComputer<FluidSystem,
1910 CartesianIndexMapper>::
1911applyNumericalAquifers_(
const GridView& gridView,
1912 const NumericalAquifers& aquifer,
1913 const bool co2store_or_h2store)
1915 const auto num_aqu_cells = aquifer.allAquiferCells();
1916 if (num_aqu_cells.empty())
return;
1919 bool oil_as_brine = co2store_or_h2store && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1920 const auto watPos = oil_as_brine? FluidSystem::oilPhaseIdx : FluidSystem::waterPhaseIdx;
1921 if (!FluidSystem::phaseIsActive(watPos)){
1922 throw std::logic_error {
"Water phase has to be active for numerical aquifer case" };
1925 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1926 auto elemIt = gridView.template begin<0>();
1927 const auto& elemEndIt = gridView.template end<0>();
1928 const auto oilPos = FluidSystem::oilPhaseIdx;
1929 const auto gasPos = FluidSystem::gasPhaseIdx;
1930 for (; elemIt != elemEndIt; ++elemIt) {
1931 const Element& element = *elemIt;
1932 const unsigned int elemIdx = elemMapper.index(element);
1933 const auto cartIx = cartesianIndexMapper_.cartesianIndex(elemIdx);
1934 const auto search = num_aqu_cells.find(cartIx);
1935 if (search != num_aqu_cells.end()) {
1937 this->sat_[watPos][elemIdx] = 1.;
1939 if (!co2store_or_h2store && FluidSystem::phaseIsActive(oilPos)) {
1940 this->sat_[oilPos][elemIdx] = 0.;
1943 if (FluidSystem::phaseIsActive(gasPos)) {
1944 this->sat_[gasPos][elemIdx] = 0.;
1946 const auto* aqu_cell = num_aqu_cells.at(cartIx);
1947 const auto msg = fmt::format(
"FOR AQUIFER CELL AT ({}, {}, {}) OF NUMERICAL "
1948 "AQUIFER {}, WATER SATURATION IS SET TO BE UNITY",
1949 aqu_cell->I+1, aqu_cell->J+1, aqu_cell->K+1, aqu_cell->aquifer_id);
1954 if (aqu_cell->init_pressure) {
1955 const Scalar pres = *(aqu_cell->init_pressure);
1956 this->pp_[watPos][elemIdx] = pres;
1957 if (FluidSystem::phaseIsActive(gasPos)) {
1958 this->pp_[gasPos][elemIdx] = pres;
1960 if (FluidSystem::phaseIsActive(oilPos)) {
1961 this->pp_[oilPos][elemIdx] = pres;
1968template<
class FluidSystem,
1971 class ElementMapper,
1972 class CartesianIndexMapper>
1974void InitialStateComputer<FluidSystem,
1978 CartesianIndexMapper>::
1979setRegionPvtIdx(
const EclipseState& eclState,
const RMap& reg)
1981 const auto& pvtnumData = eclState.fieldProps().get_int(
"PVTNUM");
1983 for (
const auto& r : reg.activeRegions()) {
1984 const auto& cells = reg.cells(r);
1985 regionPvtIdx_[r] = pvtnumData[*cells.begin()] - 1;
1989template<
class FluidSystem,
1992 class ElementMapper,
1993 class CartesianIndexMapper>
1994template<
class RMap,
class MaterialLawManager,
class Comm>
1995void InitialStateComputer<FluidSystem,
1999 CartesianIndexMapper>::
2000calcPressSatRsRv(
const RMap& reg,
2001 const std::vector<EquilRecord>& rec,
2002 MaterialLawManager& materialLawManager,
2003 const GridView& gridView,
2007 using PhaseSat = Details::PhaseSaturations<
2008 MaterialLawManager, FluidSystem, EquilReg<Scalar>,
typename RMap::CellId
2011 auto ptable = Details::PressureTable<FluidSystem, EquilReg<Scalar>>{ grav, this->num_pressure_points_ };
2012 auto psat = PhaseSat { materialLawManager, this->swatInit_ };
2013 auto vspan = std::array<Scalar, 2>{};
2015 std::vector<int> regionIsEmpty(rec.size(), 0);
2016 for (std::size_t r = 0; r < rec.size(); ++r) {
2017 const auto& cells = reg.cells(r);
2021 const auto acc = rec[r].initializationTargetAccuracy();
2025 if (cells.empty()) {
2026 regionIsEmpty[r] = 1;
2029 const auto eqreg = EquilReg {
2030 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
2031 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2034 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2035 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2036 ptable.equilibrate(eqreg, vspan);
2039 this->equilibrateTiltedFaultBlock(cells, eqreg, gridView, acc, ptable, psat);
2041 else if (acc == 0) {
2042 if (cells.empty()) {
2043 regionIsEmpty[r] = 1;
2046 const auto eqreg = EquilReg {
2047 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
2048 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2050 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2051 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2052 ptable.equilibrate(eqreg, vspan);
2054 this->equilibrateCellCentres(cells, eqreg, ptable, psat);
2057 if (cells.empty()) {
2058 regionIsEmpty[r] = 1;
2061 const auto eqreg = EquilReg {
2062 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
2063 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2065 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2066 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2067 ptable.equilibrate(eqreg, vspan);
2069 this->equilibrateHorizontal(cells, eqreg, -acc, ptable, psat);
2072 comm.min(regionIsEmpty.data(),regionIsEmpty.size());
2073 if (comm.rank() == 0) {
2074 for (std::size_t r = 0; r < rec.size(); ++r) {
2075 if (regionIsEmpty[r])
2077 +
" has no active cells");
2082template<
class FluidSystem,
2085 class ElementMapper,
2086 class CartesianIndexMapper>
2087template<
class CellRange,
class EquilibrationMethod>
2088void InitialStateComputer<FluidSystem,
2092 CartesianIndexMapper>::
2093cellLoop(
const CellRange& cells,
2094 EquilibrationMethod&& eqmethod)
2096 const auto oilPos = FluidSystem::oilPhaseIdx;
2097 const auto gasPos = FluidSystem::gasPhaseIdx;
2098 const auto watPos = FluidSystem::waterPhaseIdx;
2100 const auto oilActive = FluidSystem::phaseIsActive(oilPos);
2101 const auto gasActive = FluidSystem::phaseIsActive(gasPos);
2102 const auto watActive = FluidSystem::phaseIsActive(watPos);
2104 auto pressures = Details::PhaseQuantityValue<Scalar>{};
2105 auto saturations = Details::PhaseQuantityValue<Scalar>{};
2110 for (
const auto& cell : cells) {
2111 eqmethod(cell, pressures, saturations, Rs, Rv, Rvw);
2114 this->pp_ [oilPos][cell] = pressures.oil;
2115 this->sat_[oilPos][cell] = saturations.oil;
2119 this->pp_ [gasPos][cell] = pressures.gas;
2120 this->sat_[gasPos][cell] = saturations.gas;
2124 this->pp_ [watPos][cell] = pressures.water;
2125 this->sat_[watPos][cell] = saturations.water;
2128 if (oilActive && gasActive) {
2129 this->rs_[cell] =
Rs;
2130 this->rv_[cell] =
Rv;
2133 if (watActive && gasActive) {
2134 this->rvw_[cell] =
Rvw;
2139template<
class FluidSystem,
2142 class ElementMapper,
2143 class CartesianIndexMapper>
2144template<
class CellRange,
class PressTable,
class PhaseSat>
2145void InitialStateComputer<FluidSystem,
2149 CartesianIndexMapper>::
2150equilibrateCellCentres(
const CellRange& cells,
2151 const EquilReg<Scalar>& eqreg,
2152 const PressTable& ptable,
2155 using CellPos =
typename PhaseSat::Position;
2156 using CellID = std::remove_cv_t<std::remove_reference_t<
2157 decltype(std::declval<CellPos>().cell)>>;
2158 this->cellLoop(cells, [
this, &eqreg, &ptable, &psat]
2160 Details::PhaseQuantityValue<Scalar>& pressures,
2161 Details::PhaseQuantityValue<Scalar>& saturations,
2164 Scalar& Rvw) ->
void
2166 const auto pos = CellPos {
2167 cell, cellCenterDepth_[cell]
2170 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2171 pressures = psat.correctedPhasePressures();
2173 const auto temp = this->temperature_[cell];
2175 Rs = eqreg.dissolutionCalculator()
2176 (pos.depth, pressures.oil, temp, saturations.gas);
2178 Rv = eqreg.evaporationCalculator()
2179 (pos.depth, pressures.gas, temp, saturations.oil);
2181 Rvw = eqreg.waterEvaporationCalculator()
2182 (pos.depth, pressures.gas, temp, saturations.water);
2186template<
class FluidSystem,
2189 class ElementMapper,
2190 class CartesianIndexMapper>
2191template<
class CellRange,
class PressTable,
class PhaseSat>
2192void InitialStateComputer<FluidSystem,
2196 CartesianIndexMapper>::
2197equilibrateHorizontal(
const CellRange& cells,
2198 const EquilReg<Scalar>& eqreg,
2200 const PressTable& ptable,
2203 using CellPos =
typename PhaseSat::Position;
2204 using CellID = std::remove_cv_t<std::remove_reference_t<
2205 decltype(std::declval<CellPos>().cell)>>;
2207 this->cellLoop(cells, [
this, acc, &eqreg, &ptable, &psat]
2209 Details::PhaseQuantityValue<Scalar>& pressures,
2210 Details::PhaseQuantityValue<Scalar>& saturations,
2213 Scalar& Rvw) ->
void
2216 saturations.reset();
2218 Scalar totfrac = 0.0;
2220 const auto pos = CellPos { cell, depth };
2222 saturations.axpy(psat.deriveSaturations(pos, eqreg, ptable), frac);
2223 pressures .axpy(psat.correctedPhasePressures(), frac);
2229 saturations /= totfrac;
2230 pressures /= totfrac;
2233 const auto pos = CellPos {
2234 cell, cellCenterDepth_[cell]
2237 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2238 pressures = psat.correctedPhasePressures();
2241 const auto temp = this->temperature_[cell];
2242 const auto cz = cellCenterDepth_[cell];
2244 Rs = eqreg.dissolutionCalculator()
2245 (cz, pressures.oil, temp, saturations.gas);
2247 Rv = eqreg.evaporationCalculator()
2248 (cz, pressures.gas, temp, saturations.oil);
2250 Rvw = eqreg.waterEvaporationCalculator()
2251 (cz, pressures.gas, temp, saturations.water);
2255template<
class Flu
idSystem,
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper>
2256template<
class CellRange,
class PressTable,
class PhaseSat>
2257void InitialStateComputer<FluidSystem, Grid, GridView, ElementMapper, CartesianIndexMapper>::
2258equilibrateTiltedFaultBlockSimple(
const CellRange& cells,
2259 const EquilReg<Scalar>& eqreg,
2260 const GridView& gridView,
2262 const PressTable& ptable,
2265 using CellPos =
typename PhaseSat::Position;
2266 using CellID = std::remove_cv_t<std::remove_reference_t<
2267 decltype(std::declval<CellPos>().cell)>>;
2269 this->cellLoop(cells, [
this, acc, &eqreg, &ptable, &psat, &gridView]
2271 Details::PhaseQuantityValue<Scalar>& pressures,
2272 Details::PhaseQuantityValue<Scalar>& saturations,
2275 Scalar& Rvw) ->
void
2278 saturations.reset();
2279 Scalar totalWeight = 0.0;
2282 const auto& [zmin, zmax] = cellZMinMax_[cell];
2283 const Scalar cellThickness = zmax - zmin;
2284 const Scalar halfThickness = cellThickness / 2.0;
2287 Scalar dipAngle, dipAzimuth;
2291 std::array<Scalar, 3> referencePoint = {
2292 cellCenterXY_[cell].first,
2293 cellCenterXY_[cell].second,
2294 cellCenterDepth_[cell]
2298 const int numLevelsPerHalf = std::min(20, acc);
2301 std::vector<std::pair<Scalar, Scalar>> levels;
2304 for (
int side = 0; side < 2; ++side) {
2305 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2307 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2309 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2313 Scalar crossSectionWeight = (halfThickness / numLevelsPerHalf);
2316 if (std::abs(dipAngle) > 1e-10) {
2317 crossSectionWeight /= std::cos(dipAngle);
2320 levels.emplace_back(depth, crossSectionWeight);
2324 for (
const auto& [depth, weight] : levels) {
2326 const auto& [x, y] = cellCenterXY_[cell];
2328 depth, x, y, dipAngle, dipAzimuth, referencePoint);
2330 const auto pos = CellPos{cell, tvd};
2332 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2333 auto localPressures = psat.correctedPhasePressures();
2336 saturations.axpy(localSaturations, weight);
2337 pressures.axpy(localPressures, weight);
2338 totalWeight += weight;
2342 if (totalWeight > 1e-10) {
2343 saturations /= totalWeight;
2344 pressures /= totalWeight;
2347 const auto& [x, y] = cellCenterXY_[cell];
2349 cellCenterDepth_[cell], x, y, dipAngle, dipAzimuth, referencePoint);
2350 const auto pos = CellPos{cell, tvdCenter};
2351 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2352 pressures = psat.correctedPhasePressures();
2356 const auto temp = this->temperature_[cell];
2357 const auto& [x, y] = cellCenterXY_[cell];
2359 cellCenterDepth_[cell], x, y, dipAngle, dipAzimuth, referencePoint);
2361 Rs = eqreg.dissolutionCalculator()(tvdCenter, pressures.oil, temp, saturations.gas);
2362 Rv = eqreg.evaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.oil);
2363 Rvw = eqreg.waterEvaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.water);
2367template<
class Flu
idSystem,
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper>
2368template<
class CellRange,
class PressTable,
class PhaseSat>
2369void InitialStateComputer<FluidSystem, Grid, GridView, ElementMapper, CartesianIndexMapper>::
2370equilibrateTiltedFaultBlock(
const CellRange& cells,
2371 const EquilReg<Scalar>& eqreg,
2372 const GridView& gridView,
2374 const PressTable& ptable,
2377 using CellPos =
typename PhaseSat::Position;
2378 using CellID = std::remove_cv_t<std::remove_reference_t<
2379 decltype(std::declval<CellPos>().cell)>>;
2381 std::vector<typename GridView::template Codim<0>::Entity> entityMap(gridView.size(0));
2382 for (
const auto& entity : entities(gridView, Dune::Codim<0>())) {
2383 CellID idx = gridView.indexSet().index(entity);
2384 entityMap[idx] = entity;
2388 auto polygonArea = [](
const std::vector<std::array<Scalar, 2>>& pts) {
2389 if (pts.size() < 3)
return Scalar(0);
2391 for (
size_t i = 0; i < pts.size(); ++i) {
2392 size_t j = (i + 1) % pts.size();
2393 area += pts[i][0] * pts[j][1] - pts[j][0] * pts[i][1];
2395 return std::abs(area) * Scalar(0.5);
2399 auto computeCrossSectionArea = [&](
const CellID cell, Scalar depth) -> Scalar {
2401 const auto& entity = entityMap[cell];
2402 const auto& geometry = entity.geometry();
2403 const int numCorners = geometry.corners();
2405 std::vector<std::array<Scalar, 3>> corners(numCorners);
2406 for (
int i = 0; i < numCorners; ++i) {
2407 const auto& corner = geometry.corner(i);
2408 corners[i] = {
static_cast<Scalar
>(corner[0]),
static_cast<Scalar
>(corner[1]),
static_cast<Scalar
>(corner[2])};
2412 std::vector<std::array<Scalar, 2>> intersectionPoints;
2413 const Scalar tol = 1e-10;
2416 for (
size_t i = 0; i < corners.size(); ++i) {
2417 for (
size_t j = i + 1; j < corners.size(); ++j) {
2418 Scalar za = corners[i][2];
2419 Scalar zb = corners[j][2];
2421 if ((za - depth) * (zb - depth) <= 0.0 && std::abs(za - zb) > tol) {
2423 Scalar t = (depth - za) / (zb - za);
2424 Scalar x = corners[i][0] + t * (corners[j][0] - corners[i][0]);
2425 Scalar y = corners[i][1] + t * (corners[j][1] - corners[i][1]);
2426 intersectionPoints.push_back({x, y});
2432 if (intersectionPoints.size() > 1) {
2433 auto pointsEqual = [tol](
const std::array<Scalar, 2>& a,
const std::array<Scalar, 2>& b) {
2434 return std::abs(a[0] - b[0]) < tol && std::abs(a[1] - b[1]) < tol;
2437 intersectionPoints.erase(
2438 std::unique(intersectionPoints.begin(), intersectionPoints.end(), pointsEqual),
2439 intersectionPoints.end()
2443 if (intersectionPoints.size() < 3) {
2449 Scalar cx = 0, cy = 0;
2450 for (
const auto& p : intersectionPoints) {
2451 cx += p[0]; cy += p[1];
2453 cx /= intersectionPoints.size();
2454 cy /= intersectionPoints.size();
2457 auto angleCompare = [cx, cy](
const std::array<Scalar, 2>& a,
const std::array<Scalar, 2>& b) {
2458 return std::atan2(a[1] - cy, a[0] - cx) < std::atan2(b[1] - cy, b[0] - cx);
2461 std::ranges::sort(intersectionPoints, angleCompare);
2463 return polygonArea(intersectionPoints);
2465 }
catch (
const std::exception& e) {
2470 auto cellProcessor = [
this, acc, &eqreg, &ptable, &psat, &computeCrossSectionArea]
2472 Details::PhaseQuantityValue<Scalar>& pressures,
2473 Details::PhaseQuantityValue<Scalar>& saturations,
2476 Scalar&
Rvw) ->
void
2479 saturations.reset();
2480 Scalar totalWeight = 0.0;
2482 const auto& zmin = this->cellZMinMax_[cell].first;
2483 const auto& zmax = this->cellZMinMax_[cell].second;
2484 const Scalar cellThickness = zmax - zmin;
2485 const Scalar halfThickness = cellThickness / 2.0;
2488 Scalar dipAngle, dipAzimuth;
2492 std::array<Scalar, 3> referencePoint = {
2493 this->cellCenterXY_[cell].first,
2494 this->cellCenterXY_[cell].second,
2495 cellCenterDepth_[cell]
2499 const int numLevelsPerHalf = std::min(20, acc);
2502 std::vector<std::pair<Scalar, Scalar>> levels;
2505 for (
int side = 0; side < 2; ++side) {
2506 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2508 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2510 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2513 Scalar crossSectionArea = computeCrossSectionArea(cell, depth);
2516 Scalar weight = crossSectionArea * (halfThickness / numLevelsPerHalf);
2518 levels.emplace_back(depth, weight);
2523 bool hasValidAreas =
false;
2524 for (
const auto& level : levels) {
2525 if (level.second > 1e-10) {
2526 hasValidAreas =
true;
2531 if (!hasValidAreas) {
2534 for (
int side = 0; side < 2; ++side) {
2535 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2536 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2537 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2538 Scalar weight = (halfThickness / numLevelsPerHalf);
2539 if (std::abs(dipAngle) > 1e-10) {
2540 weight /= std::cos(dipAngle);
2542 levels.emplace_back(depth, weight);
2547 for (
const auto& level : levels) {
2548 Scalar depth = level.first;
2549 Scalar weight = level.second;
2552 const auto& xy = this->cellCenterXY_[cell];
2554 depth, xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2556 const auto pos = CellPos{cell, tvd};
2558 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2559 auto localPressures = psat.correctedPhasePressures();
2562 saturations.axpy(localSaturations, weight);
2563 pressures.axpy(localPressures, weight);
2564 totalWeight += weight;
2567 if (totalWeight > 1e-10) {
2568 saturations /= totalWeight;
2569 pressures /= totalWeight;
2572 const auto& xy = this->cellCenterXY_[cell];
2574 this->cellCenterDepth_[cell], xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2575 const auto pos = CellPos{cell, tvdCenter};
2576 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2577 pressures = psat.correctedPhasePressures();
2581 const auto temp = this->temperature_[cell];
2582 const auto& xy = this->cellCenterXY_[cell];
2584 this->cellCenterDepth_[cell], xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2586 Rs = eqreg.dissolutionCalculator()(tvdCenter, pressures.oil, temp, saturations.gas);
2587 Rv = eqreg.evaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.oil);
2588 Rvw = eqreg.waterEvaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.water);
2591 this->cellLoop(cells, cellProcessor);
#define OPM_END_PARALLEL_TRY_CATCH(prefix, comm)
Catch exception and throw in a parallel try-catch clause.
Definition: DeferredLoggingErrorHelpers.hpp:192
#define OPM_BEGIN_PARALLEL_TRY_CATCH()
Macro to setup the try of a parallel try-catch.
Definition: DeferredLoggingErrorHelpers.hpp:158
Auxiliary routines that to solve the ODEs that emerge from the hydrostatic equilibrium problem.
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:325
Routines that actually solve the ODEs that emerge from the hydrostatic equilibrium problem.
Definition: InitStateEquil.hpp:704
Definition: InitStateEquil.hpp:151
Gas(const TabulatedFunction &tempVdTable, const RV &rv, const RVW &rvw, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:476
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:492
Definition: InitStateEquil.hpp:126
Oil(const TabulatedFunction &tempVdTable, const RS &rs, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:428
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:442
Definition: InitStateEquil.hpp:101
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:402
Water(const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:388
Definition: InitStateEquil.hpp:399
const PhaseQuantityValue< Scalar > & deriveSaturations(const Position &x, const Region ®, const PTable &ptable)
Definition: InitStateEquil_impl.hpp:734
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< Scalar > &swatInit)
Definition: InitStateEquil_impl.hpp:710
Definition: InitStateEquil.hpp:180
PressureTable & operator=(const PressureTable &rhs)
Definition: InitStateEquil_impl.hpp:1155
Scalar water(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1235
Scalar gas(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1224
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition: InitStateEquil_impl.hpp:1206
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition: InitStateEquil_impl.hpp:1199
Scalar oil(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1214
std::array< Scalar, 2 > VSpan
Definition: InitStateEquil.hpp:183
bool oilActive() const
Predicate for whether or not oil is an active phase.
Definition: InitStateEquil_impl.hpp:1192
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:182
void equilibrate(const Region ®, const VSpan &span)
Definition: InitStateEquil_impl.hpp:1181
PressureTable(const Scalar gravity, const int samplePoints=2000)
Definition: InitStateEquil_impl.hpp:1125
Definition: InitStateEquil.hpp:80
Scalar operator()(const Scalar x) const
Definition: InitStateEquil_impl.hpp:354
RK4IVP(const RHS &f, const std::array< Scalar, 2 > &span, const Scalar y0, const int N)
Definition: InitStateEquil_impl.hpp:319
Definition: EquilibrationHelpers.hpp:135
Definition: EquilibrationHelpers.hpp:216
Definition: EquilibrationHelpers.hpp:269
Definition: EquilibrationHelpers.hpp:612
Definition: EquilibrationHelpers.hpp:162
Definition: EquilibrationHelpers.hpp:322
Definition: EquilibrationHelpers.hpp:376
std::vector< EquilRecord > getEquil(const EclipseState &state)
Definition: InitStateEquil_impl.hpp:1429
std::vector< int > equilnum(const EclipseState &eclipseState, const GridView &gridview)
Definition: InitStateEquil_impl.hpp:1443
std::pair< Scalar, Scalar > cellZMinMax(const Element &element)
Definition: InitStateEquil_impl.hpp:186
Scalar cellCenterDepth(const Element &element)
Definition: InitStateEquil_impl.hpp:133
std::pair< Scalar, Scalar > cellZSpan(const Element &element)
Definition: InitStateEquil_impl.hpp:167
CellCornerData< Scalar > getCellCornerXY(const Element &element)
Definition: InitStateEquil_impl.hpp:256
void verticalExtent(const CellRange &cells, const std::vector< std::pair< Scalar, Scalar > > &cellZMinMax, const Parallel::Communication &comm, std::array< Scalar, 2 > &span)
Definition: InitStateEquil_impl.hpp:69
std::pair< Scalar, Scalar > cellCenterXY(const Element &element)
Definition: InitStateEquil_impl.hpp:148
Scalar calculateTrueVerticalDepth(Scalar z, Scalar x, Scalar y, Scalar dipAngle, Scalar dipAzimuth, const std::array< Scalar, 3 > &referencePoint)
Definition: InitStateEquil_impl.hpp:280
void subdivisionCentrePoints(const Scalar left, const Scalar right, const int numIntervals, std::vector< std::pair< Scalar, Scalar > > &subdiv)
Definition: InitStateEquil_impl.hpp:94
std::vector< std::pair< Scalar, Scalar > > horizontalSubdivision(const CellID cell, const std::pair< Scalar, Scalar > topbot, const int numIntervals)
Definition: InitStateEquil_impl.hpp:112
void computeBlockDip(const CellCornerData< Scalar > &cellCorners, Scalar &dipAngle, Scalar &dipAzimuth)
Definition: InitStateEquil_impl.hpp:205
Dune::Communication< MPIComm > Communication
Definition: ParallelCommunication.hpp:30
Definition: blackoilbioeffectsmodules.hh:45
std::string to_string(const ConvergenceReport::ReservoirFailure::Type t)
Definition: InitStateEquil.hpp:61
std::array< Scalar, 8 > X
Definition: InitStateEquil.hpp:62
std::array< Scalar, 8 > Y
Definition: InitStateEquil.hpp:63
std::array< Scalar, 8 > Z
Definition: InitStateEquil.hpp:64
Simple set of per-phase (named by primary component) quantities.
Definition: InitStateEquil.hpp:351
Definition: InitStateEquil.hpp:405