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/RsvdTable.hpp>
34#include <opm/input/eclipse/EclipseState/Tables/RvvdTable.hpp>
35#include <opm/input/eclipse/EclipseState/Tables/RvwvdTable.hpp>
36#include <opm/input/eclipse/EclipseState/Tables/PbvdTable.hpp>
37#include <opm/input/eclipse/EclipseState/Tables/PdvdTable.hpp>
38#include <opm/input/eclipse/EclipseState/Tables/SaltvdTable.hpp>
39#include <opm/input/eclipse/EclipseState/Tables/RtempvdTable.hpp>
41#include <opm/input/eclipse/EclipseState/Tables/SaltpvdTable.hpp>
43#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
44#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
50#include <fmt/format.h>
64template <
typename CellRange,
class Scalar>
66 const std::vector<std::pair<Scalar, Scalar>>&
cellZMinMax,
68 std::array<Scalar,2>& span)
70 span[0] = std::numeric_limits<Scalar>::max();
71 span[1] = std::numeric_limits<Scalar>::lowest();
81 for (
const auto& cell : cells) {
85 span[0] = comm.min(span[0]);
86 span[1] = comm.max(span[1]);
92 const int numIntervals,
93 std::vector<std::pair<Scalar, Scalar>>& subdiv)
95 const auto h = (right - left) / numIntervals;
98 for (
auto i = 0*numIntervals; i < numIntervals; ++i) {
99 const auto start = end;
100 end = left + (i + 1)*h;
102 subdiv.emplace_back((start + end) / 2, h);
106template <
typename CellID,
typename Scalar>
107std::vector<std::pair<Scalar, Scalar>>
109 const std::pair<Scalar, Scalar> topbot,
110 const int numIntervals)
112 auto subdiv = std::vector<std::pair<Scalar, Scalar>>{};
113 subdiv.reserve(2 * numIntervals);
115 if (topbot.first > topbot.second) {
116 throw std::out_of_range {
117 "Negative thickness (inverted top/bottom faces) in cell "
123 2*numIntervals, subdiv);
128template <
class Scalar,
class Element>
131 typedef typename Element::Geometry Geometry;
132 static constexpr int zCoord = Element::dimension - 1;
135 const Geometry& geometry = element.geometry();
136 const int corners = geometry.corners();
137 for (
int i=0; i < corners; ++i)
138 zz += geometry.corner(i)[zCoord];
143template <
class Scalar,
class Element>
146 typedef typename Element::Geometry Geometry;
147 static constexpr int xCoord = Element::dimension - 3;
148 static constexpr int yCoord = Element::dimension - 2;
153 const Geometry& geometry = element.geometry();
154 const int corners = geometry.corners();
155 for (
int i=0; i < corners; ++i) {
156 xx += geometry.corner(i)[xCoord];
157 yy += geometry.corner(i)[yCoord];
159 return std::make_pair(xx/corners, yy/corners);
162template <
class Scalar,
class Element>
163std::pair<Scalar,Scalar>
cellZSpan(
const Element& element)
165 typedef typename Element::Geometry Geometry;
166 static constexpr int zCoord = Element::dimension - 1;
170 const Geometry& geometry = element.geometry();
171 const int corners = geometry.corners();
172 assert(corners == 8);
173 for (
int i=0; i < 4; ++i)
174 bot += geometry.corner(i)[zCoord];
175 for (
int i=4; i < corners; ++i)
176 top += geometry.corner(i)[zCoord];
178 return std::make_pair(bot/4, top/4);
181template <
class Scalar,
class Element>
184 typedef typename Element::Geometry Geometry;
185 static constexpr int zCoord = Element::dimension - 1;
186 const Geometry& geometry = element.geometry();
187 const int corners = geometry.corners();
188 assert(corners == 8);
189 auto min = std::numeric_limits<Scalar>::max();
190 auto max = std::numeric_limits<Scalar>::lowest();
193 for (
int i=0; i < corners; ++i) {
194 min = std::min(min,
static_cast<Scalar
>(geometry.corner(i)[zCoord]));
195 max = std::max(max,
static_cast<Scalar
>(geometry.corner(i)[zCoord]));
197 return std::make_pair(min, max);
200template<
class Scalar>
202 Scalar& dipAngle, Scalar& dipAzimuth)
204 const auto& Xc = cellCorners.
X;
205 const auto& Yc = cellCorners.
Y;
206 const auto& Zc = cellCorners.
Z;
208 Scalar v1x = Xc[1] - Xc[0];
209 Scalar v1y = Yc[1] - Yc[0];
210 Scalar v1z = Zc[1] - Zc[0];
212 Scalar v2x = Xc[2] - Xc[0];
213 Scalar v2y = Yc[2] - Yc[0];
214 Scalar v2z = Zc[2] - Zc[0];
217 Scalar nx = v1y * v2z - v1z * v2y;
218 Scalar ny = v1z * v2x - v1x * v2z;
219 Scalar nz = v1x * v2y - v1y * v2x;
222 Scalar norm = std::hypot(nx, ny, nz);
230 dipAngle = std::acos(std::abs(nz));
233 if (std::abs(nx) > 1e-10 || std::abs(ny) > 1e-10) {
234 dipAzimuth = std::atan2(ny, nx);
236 dipAzimuth = std::fmod(dipAzimuth + 2*M_PI, 2*M_PI);
242 const Scalar maxDip = M_PI/2 - 1e-6;
243 dipAngle = std::min(dipAngle, maxDip);
251template <
class Scalar,
class Element>
254 typedef typename Element::Geometry Geometry;
255 const Geometry& geometry = element.geometry();
256 static constexpr int zCoord = Element::dimension - 1;
257 static constexpr int yCoord = Element::dimension - 2;
258 static constexpr int xCoord = Element::dimension - 3;
259 const int corners = geometry.corners();
260 assert(corners == 8);
261 std::array<Scalar, 8> X {};
262 std::array<Scalar, 8> Y {};
263 std::array<Scalar, 8> Z {};
265 for (
int i = 0; i < corners; ++i) {
266 auto corner = geometry.corner(i);
267 X[i] = corner[xCoord];
268 Y[i] = corner[yCoord];
269 Z[i] = corner[zCoord];
275template<
class Scalar>
277 Scalar dipAngle, Scalar dipAzimuth,
278 const std::array<Scalar, 3>& referencePoint)
285 Scalar dx = x - referencePoint[0];
286 Scalar dy = y - referencePoint[1];
287 Scalar dz = z - referencePoint[2];
290 if (std::abs(dipAngle) < 1e-10) {
291 return referencePoint[2] + dz;
295 Scalar pointAzimuth = std::atan2(dy, dx);
298 Scalar azimuthDiff = pointAzimuth - dipAzimuth;
301 Scalar lateralDist = std::hypot(dx, dy);
304 Scalar lateralInDipDir = lateralDist * std::cos(azimuthDiff);
309 Scalar tvd = referencePoint[2] + dz * std::cos(dipAngle) + lateralInDipDir * std::sin(dipAngle);
314template<
class Scalar,
class RHS>
316 const std::array<Scalar,2>& span,
322 const Scalar h = stepsize();
323 const Scalar h2 = h / 2;
324 const Scalar h6 = h / 6;
330 f_.push_back(f(span_[0], y0));
332 for (
int i = 0; i < N; ++i) {
333 const Scalar x = span_[0] + i*h;
334 const Scalar y = y_.back();
336 const Scalar k1 = f_[i];
337 const Scalar k2 = f(x + h2, y + h2*k1);
338 const Scalar k3 = f(x + h2, y + h2*k2);
339 const Scalar k4 = f(x + h, y + h*k3);
341 y_.push_back(y + h6*(k1 + 2*(k2 + k3) + k4));
342 f_.push_back(f(x + h, y_.back()));
345 assert (y_.size() ==
typename std::vector<Scalar>::size_type(N + 1));
348template<
class Scalar,
class RHS>
354 const Scalar h = stepsize();
355 int i = (x - span_[0]) / h;
356 const Scalar t = (x - (span_[0] + i*h)) / h;
359 if (i < 0) { i = 0; }
360 if (N_ <= i) { i = N_ - 1; }
362 const Scalar y0 = y_[i], y1 = y_[i + 1];
363 const Scalar f0 = f_[i], f1 = f_[i + 1];
365 Scalar u = (1 - 2*t) * (y1 - y0);
366 u += h * ((t - 1)*f0 + t*f1);
368 u += (1 - t)*y0 + t*y1;
373template<
class Scalar,
class RHS>
377 return (span_[1] - span_[0]) / N_;
380namespace PhasePressODE {
382template<
class Flu
idSystem>
384Water(
const TabulatedFunction& tempVdTable,
385 const TabulatedFunction& saltVdTable,
386 const int pvtRegionIdx,
387 const Scalar normGrav)
388 : tempVdTable_(tempVdTable)
389 , saltVdTable_(saltVdTable)
390 , pvtRegionIdx_(pvtRegionIdx)
395template<
class Flu
idSystem>
396typename Water<FluidSystem>::Scalar
399 const Scalar press)
const
401 return this->density(depth, press) * g_;
404template<
class Flu
idSystem>
405typename Water<FluidSystem>::Scalar
408 const Scalar press)
const
411 Scalar saltConcentration = saltVdTable_.eval(depth,
true);
412 Scalar temp = tempVdTable_.eval(depth,
true);
413 Scalar rho = FluidSystem::waterPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
418 rho *= FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
422template<
class Flu
idSystem,
class RS>
424Oil(
const TabulatedFunction& tempVdTable,
426 const int pvtRegionIdx,
427 const Scalar normGrav)
428 : tempVdTable_(tempVdTable)
430 , pvtRegionIdx_(pvtRegionIdx)
435template<
class Flu
idSystem,
class RS>
436typename Oil<FluidSystem,RS>::Scalar
439 const Scalar press)
const
441 return this->density(depth, press) * g_;
444template<
class Flu
idSystem,
class RS>
445typename Oil<FluidSystem,RS>::Scalar
448 const Scalar press)
const
450 const Scalar temp = tempVdTable_.eval(depth,
true);
452 if (FluidSystem::enableDissolvedGas())
453 rs = rs_(depth, press, temp);
456 if (rs >= FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press)) {
457 bOil = FluidSystem::oilPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
460 bOil = FluidSystem::oilPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rs);
462 Scalar rho = bOil * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
463 if (FluidSystem::enableDissolvedGas()) {
464 rho += rs * bOil * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
470template<
class Flu
idSystem,
class RV,
class RVW>
472Gas(
const TabulatedFunction& tempVdTable,
475 const int pvtRegionIdx,
476 const Scalar normGrav)
477 : tempVdTable_(tempVdTable)
480 , pvtRegionIdx_(pvtRegionIdx)
485template<
class Flu
idSystem,
class RV,
class RVW>
486typename Gas<FluidSystem,RV,RVW>::Scalar
489 const Scalar press)
const
491 return this->density(depth, press) * g_;
494template<
class Flu
idSystem,
class RV,
class RVW>
495typename Gas<FluidSystem,RV,RVW>::Scalar
498 const Scalar press)
const
500 const Scalar temp = tempVdTable_.eval(depth,
true);
502 if (FluidSystem::enableVaporizedOil())
503 rv = rv_(depth, press, temp);
506 if (FluidSystem::enableVaporizedWater())
507 rvw = rvw_(depth, press, temp);
511 if (FluidSystem::enableVaporizedOil() && FluidSystem::enableVaporizedWater()) {
512 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)
513 && rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press))
515 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
517 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rv, rvw);
519 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
520 rho += rv * bGas * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_)
521 + rvw * bGas * FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
525 if (FluidSystem::enableVaporizedOil()){
526 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)) {
527 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
529 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
535 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
536 rho += rv * bGas * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
540 if (FluidSystem::enableVaporizedWater()){
541 if (rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press)) {
542 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
545 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
551 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
552 rho += rvw * bGas * FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
557 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp,
561 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
568template<
class Flu
idSystem,
class Region>
570PressureTable<FluidSystem,Region>::
571PressureFunction<ODE>::PressureFunction(
const ODE& ode,
577 this->value_[Direction::Up] = std::make_unique<Distribution>
578 (ode, VSpan {{ ic.depth, span[0] }}, ic.pressure, nsample);
580 this->value_[Direction::Down] = std::make_unique<Distribution>
581 (ode, VSpan {{ ic.depth, span[1] }}, ic.pressure, nsample);
584template<
class Flu
idSystem,
class Region>
586PressureTable<FluidSystem,Region>::
587PressureFunction<ODE>::PressureFunction(
const PressureFunction& rhs)
588 : initial_(rhs.initial_)
590 this->value_[Direction::Up] =
591 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
593 this->value_[Direction::Down] =
594 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
597template<
class Flu
idSystem,
class Region>
599typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
604 this->initial_ = rhs.initial_;
606 this->value_[Direction::Up] =
607 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
609 this->value_[Direction::Down] =
610 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
615template<
class Flu
idSystem,
class Region>
617typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
622 this->initial_ = rhs.initial_;
623 this->value_ = std::move(rhs.value_);
628template<
class Flu
idSystem,
class Region>
631PressureTable<FluidSystem,Region>::
632PressureFunction<ODE>::
633value(
const Scalar depth)
const
635 if (depth < this->initial_.depth) {
637 return (*this->value_[Direction::Up])(depth);
639 else if (depth > this->initial_.depth) {
641 return (*this->value_[Direction::Down])(depth);
645 return this->initial_.pressure;
650template<
class Flu
idSystem,
class Region>
651template<
typename PressFunc>
652void PressureTable<FluidSystem,Region>::
653checkPtr(
const PressFunc* phasePress,
654 const std::string& phaseName)
const
656 if (phasePress !=
nullptr) {
return; }
658 throw std::invalid_argument {
659 "Phase pressure function for \"" + phaseName
660 +
"\" most not be null"
664template<
class Flu
idSystem,
class Region>
665typename PressureTable<FluidSystem,Region>::Strategy
666PressureTable<FluidSystem,Region>::
667selectEquilibrationStrategy(
const Region& reg)
const
669 if (!this->oilActive()) {
670 if (reg.datum() > reg.zwoc()) {
671 return &PressureTable::equil_WOG;
673 return &PressureTable::equil_GOW;
676 if (reg.datum() > reg.zwoc()) {
677 return &PressureTable::equil_WOG;
679 else if (reg.datum() < reg.zgoc()) {
680 return &PressureTable::equil_GOW;
683 return &PressureTable::equil_OWG;
687template<
class Flu
idSystem,
class Region>
688void PressureTable<FluidSystem,Region>::
689copyInPointers(
const PressureTable& rhs)
691 if (rhs.oil_ !=
nullptr) {
692 this->oil_ = std::make_unique<OPress>(*rhs.oil_);
695 if (rhs.gas_ !=
nullptr) {
696 this->gas_ = std::make_unique<GPress>(*rhs.gas_);
699 if (rhs.wat_ !=
nullptr) {
700 this->wat_ = std::make_unique<WPress>(*rhs.wat_);
704template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
707 const std::vector<Scalar>& swatInit)
708 : matLawMgr_(matLawMgr)
709 , swatInit_ (swatInit)
713template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
716 : matLawMgr_(rhs.matLawMgr_)
717 , swatInit_ (rhs.swatInit_)
719 , press_ (rhs.press_)
722 this->setEvaluationPoint(*rhs.evalPt_.position,
724 *rhs.evalPt_.ptable);
727template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
734 this->setEvaluationPoint(x, reg, ptable);
735 this->initializePhaseQuantities();
737 if (ptable.
gasActive()) { this->deriveGasSat(); }
739 if (ptable.
waterActive()) { this->deriveWaterSat(); }
742 if (this->isOverlappingTransition()) {
743 this->fixUnphysicalTransition();
746 if (ptable.
oilActive()) { this->deriveOilSat(); }
748 this->accountForScaledSaturations();
753template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
757 const PTable& ptable)
759 this->evalPt_.position = &x;
760 this->evalPt_.region = ®
761 this->evalPt_.ptable = &ptable;
764template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
765void PhaseSaturations<MaterialLawManager,FluidSystem,Region,CellID>::
766initializePhaseQuantities()
769 this->press_.reset();
771 const auto depth = this->evalPt_.position->depth;
772 const auto& ptable = *this->evalPt_.ptable;
774 if (ptable.oilActive()) {
775 this->press_.oil = ptable.oil(depth);
778 if (ptable.gasActive()) {
779 this->press_.gas = ptable.gas(depth);
782 if (ptable.waterActive()) {
783 this->press_.water = ptable.water(depth);
787template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
788void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveOilSat()
790 this->sat_.oil = 1.0 - this->sat_.water - this->sat_.gas;
793template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
794void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveGasSat()
796 auto& sg = this->sat_.gas;
798 const auto isIncr =
true;
799 const auto oilActive = this->evalPt_.ptable->oilActive();
801 if (this->isConstCapPress(this->gasPos())) {
805 const auto gas_contact = oilActive? this->evalPt_.region->zgoc() : this->evalPt_.region->zwoc();
806 sg = this->fromDepthTable(gas_contact,
807 this->gasPos(), isIncr);
817 const auto pw = oilActive? this->press_.oil : this->press_.water;
818 const auto pcgo = this->press_.gas - pw;
819 sg = this->invertCapPress(pcgo, this->gasPos(), isIncr);
823template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
824void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveWaterSat()
826 auto& sw = this->sat_.water;
828 const auto oilActive = this->evalPt_.ptable->oilActive();
831 sw = 1.0 - this->sat_.gas;
834 const auto isIncr =
false;
836 if (this->isConstCapPress(this->waterPos())) {
840 sw = this->fromDepthTable(this->evalPt_.region->zwoc(),
841 this->waterPos(), isIncr);
853 const auto pcow = this->press_.oil - this->press_.water;
855 if (this->swatInit_.empty()) {
856 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
859 auto [swout, newSwatInit] = this->applySwatInit(pcow);
861 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
870template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
871void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
872fixUnphysicalTransition()
874 auto& sg = this->sat_.gas;
875 auto& sw = this->sat_.water;
883 const auto pcgw = this->press_.gas - this->press_.water;
884 if (! this->swatInit_.empty()) {
888 auto [swout, newSwatInit] = this->applySwatInit(pcgw, sw);
890 const auto isIncr =
false;
891 sw = this->invertCapPress(pcgw, this->waterPos(), isIncr);
898 sw = satFromSumOfPcs<FluidSystem>
899 (this->matLawMgr_, this->waterPos(), this->gasPos(),
900 this->evalPt_.position->cell, pcgw);
903 this->fluidState_.setSaturation(this->oilPos(), 1.0 - sw - sg);
904 this->fluidState_.setSaturation(this->gasPos(), sg);
905 this->fluidState_.setSaturation(this->waterPos(), this->evalPt_
906 .ptable->waterActive() ? sw : 0.0);
909 this->computeMaterialLawCapPress();
910 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
913template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
914void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
915accountForScaledSaturations()
917 const auto gasActive = this->evalPt_.ptable->gasActive();
918 const auto watActive = this->evalPt_.ptable->waterActive();
919 const auto oilActive = this->evalPt_.ptable->oilActive();
921 auto sg = gasActive? this->sat_.gas : 0.0;
922 auto sw = watActive? this->sat_.water : 0.0;
923 auto so = oilActive? this->sat_.oil : 0.0;
925 this->fluidState_.setSaturation(this->waterPos(), sw);
926 this->fluidState_.setSaturation(this->oilPos(), so);
927 this->fluidState_.setSaturation(this->gasPos(), sg);
929 const auto& scaledDrainageInfo = this->matLawMgr_
930 .oilWaterScaledEpsInfoDrainage(this->evalPt_.position->cell);
932 const auto thresholdSat = 1.0e-6;
933 if (watActive && ((sw + thresholdSat) > scaledDrainageInfo.Swu)) {
937 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swu);
939 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swu);
940 }
else if (gasActive) {
941 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swu);
943 sw = scaledDrainageInfo.Swu;
944 this->computeMaterialLawCapPress();
948 this->press_.oil = this->press_.water + this->materialLawCapPressOilWater();
951 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
955 if (gasActive && ((sg + thresholdSat) > scaledDrainageInfo.Sgu)) {
959 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgu);
961 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgu);
962 }
else if (watActive) {
963 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgu);
965 sg = scaledDrainageInfo.Sgu;
966 this->computeMaterialLawCapPress();
970 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
973 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
977 if (watActive && ((sw - thresholdSat) < scaledDrainageInfo.Swl)) {
981 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swl);
983 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swl);
984 }
else if (gasActive) {
985 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swl);
987 sw = scaledDrainageInfo.Swl;
988 this->computeMaterialLawCapPress();
992 this->press_.water = this->press_.oil - this->materialLawCapPressOilWater();
995 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
999 if (gasActive && ((sg - thresholdSat) < scaledDrainageInfo.Sgl)) {
1003 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgl);
1005 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgl);
1006 }
else if (watActive) {
1007 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgl);
1009 sg = scaledDrainageInfo.Sgl;
1010 this->computeMaterialLawCapPress();
1014 this->press_.gas = this->press_.oil + this->materialLawCapPressGasOil();
1017 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
1022template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1023std::pair<typename FluidSystem::Scalar, bool>
1024PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1025applySwatInit(
const Scalar pcow)
1027 return this->applySwatInit(pcow, this->swatInit_[this->evalPt_.position->cell]);
1030template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1031std::pair<typename FluidSystem::Scalar, bool>
1032PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1033applySwatInit(
const Scalar pcow,
const Scalar sw)
1035 return this->matLawMgr_.applySwatinit(this->evalPt_.position->cell, pcow, sw);
1038template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1039void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1040computeMaterialLawCapPress()
1042 const auto& matParams = this->matLawMgr_
1043 .materialLawParams(this->evalPt_.position->cell);
1045 this->matLawCapPress_.fill(0.0);
1046 MaterialLaw::capillaryPressures(this->matLawCapPress_,
1047 matParams, this->fluidState_);
1050template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1051typename FluidSystem::Scalar
1052PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1053materialLawCapPressGasOil()
const
1055 return this->matLawCapPress_[this->oilPos()]
1056 + this->matLawCapPress_[this->gasPos()];
1059template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1060typename FluidSystem::Scalar
1061PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1062materialLawCapPressOilWater()
const
1064 return this->matLawCapPress_[this->oilPos()]
1065 - this->matLawCapPress_[this->waterPos()];
1068template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1069typename FluidSystem::Scalar
1070PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1071materialLawCapPressGasWater()
const
1073 return this->matLawCapPress_[this->gasPos()]
1074 - this->matLawCapPress_[this->waterPos()];
1077template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1078bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1079isConstCapPress(
const PhaseIdx phaseIdx)
const
1081 return isConstPc<FluidSystem>
1082 (this->matLawMgr_, phaseIdx, this->evalPt_.position->cell);
1085template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1086bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1087isOverlappingTransition()
const
1089 return this->evalPt_.ptable->gasActive()
1090 && this->evalPt_.ptable->waterActive()
1091 && ((this->sat_.gas + this->sat_.water) > 1.0);
1094template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1095typename FluidSystem::Scalar
1096PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1097fromDepthTable(
const Scalar contactdepth,
1098 const PhaseIdx phasePos,
1099 const bool isincr)
const
1101 return satFromDepth<FluidSystem>
1102 (this->matLawMgr_, this->evalPt_.position->depth,
1103 contactdepth,
static_cast<int>(phasePos),
1104 this->evalPt_.position->cell, isincr);
1107template <
class MaterialLawManager,
class Flu
idSystem,
class Region,
typename CellID>
1108typename FluidSystem::Scalar
1109PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1110invertCapPress(
const Scalar pc,
1111 const PhaseIdx phasePos,
1112 const bool isincr)
const
1114 return satFromPc<FluidSystem>
1115 (this->matLawMgr_,
static_cast<int>(phasePos),
1116 this->evalPt_.position->cell, pc, isincr);
1119template<
class Flu
idSystem,
class Region>
1122 const int samplePoints)
1124 , nsample_(samplePoints)
1128template <
class Flu
idSystem,
class Region>
1131 : gravity_(rhs.gravity_)
1132 , nsample_(rhs.nsample_)
1134 this->copyInPointers(rhs);
1137template <
class Flu
idSystem,
class Region>
1140 : gravity_(rhs.gravity_)
1141 , nsample_(rhs.nsample_)
1142 , oil_ (std::move(rhs.oil_))
1143 , gas_ (std::move(rhs.gas_))
1144 , wat_ (std::move(rhs.wat_))
1148template <
class Flu
idSystem,
class Region>
1153 this->gravity_ = rhs.gravity_;
1154 this->nsample_ = rhs.nsample_;
1155 this->copyInPointers(rhs);
1160template <
class Flu
idSystem,
class Region>
1165 this->gravity_ = rhs.gravity_;
1166 this->nsample_ = rhs.nsample_;
1168 this->oil_ = std::move(rhs.oil_);
1169 this->gas_ = std::move(rhs.gas_);
1170 this->wat_ = std::move(rhs.wat_);
1175template <
class Flu
idSystem,
class Region>
1181 auto equil = this->selectEquilibrationStrategy(reg);
1183 (this->*equil)(reg, span);
1186template <
class Flu
idSystem,
class Region>
1190 return FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1193template <
class Flu
idSystem,
class Region>
1197 return FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
1200template <
class Flu
idSystem,
class Region>
1204 return FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1207template <
class Flu
idSystem,
class Region>
1208typename FluidSystem::Scalar
1212 this->checkPtr(this->oil_.get(),
"OIL");
1214 return this->oil_->value(depth);
1217template <
class Flu
idSystem,
class Region>
1218typename FluidSystem::Scalar
1222 this->checkPtr(this->gas_.get(),
"GAS");
1224 return this->gas_->value(depth);
1228template <
class Flu
idSystem,
class Region>
1229typename FluidSystem::Scalar
1233 this->checkPtr(this->wat_.get(),
"WATER");
1235 return this->wat_->value(depth);
1238template <
class Flu
idSystem,
class Region>
1240equil_WOG(
const Region& reg,
const VSpan& span)
1245 if (! this->waterActive()) {
1246 throw std::invalid_argument {
1247 "Don't know how to interpret EQUIL datum depth in "
1248 "WATER zone in model without active water phase"
1253 const auto ic =
typename WPress::InitCond {
1254 reg.datum(), reg.pressure()
1257 this->makeWatPressure(ic, reg, span);
1260 if (this->oilActive()) {
1262 const auto ic =
typename OPress::InitCond {
1264 this->water(reg.zwoc()) + reg.pcowWoc()
1267 this->makeOilPressure(ic, reg, span);
1270 if (this->gasActive() && this->oilActive()) {
1272 const auto ic =
typename GPress::InitCond {
1274 this->oil(reg.zgoc()) + reg.pcgoGoc()
1277 this->makeGasPressure(ic, reg, span);
1278 }
else if (this->gasActive() && !this->oilActive()) {
1280 const auto ic =
typename GPress::InitCond {
1282 this->water(reg.zwoc()) + reg.pcowWoc()
1284 this->makeGasPressure(ic, reg, span);
1288template <
class Flu
idSystem,
class Region>
1289void PressureTable<FluidSystem, Region>::
1290equil_GOW(
const Region& reg,
const VSpan& span)
1295 if (! this->gasActive()) {
1296 throw std::invalid_argument {
1297 "Don't know how to interpret EQUIL datum depth in "
1298 "GAS zone in model without active gas phase"
1303 const auto ic =
typename GPress::InitCond {
1304 reg.datum(), reg.pressure()
1307 this->makeGasPressure(ic, reg, span);
1310 if (this->oilActive()) {
1312 const auto ic =
typename OPress::InitCond {
1314 this->gas(reg.zgoc()) - reg.pcgoGoc()
1316 this->makeOilPressure(ic, reg, span);
1319 if (this->waterActive() && this->oilActive()) {
1321 const auto ic =
typename WPress::InitCond {
1323 this->oil(reg.zwoc()) - reg.pcowWoc()
1326 this->makeWatPressure(ic, reg, span);
1327 }
else if (this->waterActive() && !this->oilActive()) {
1329 const auto ic =
typename WPress::InitCond {
1331 this->gas(reg.zwoc()) - reg.pcowWoc()
1333 this->makeWatPressure(ic, reg, span);
1337template <
class Flu
idSystem,
class Region>
1338void PressureTable<FluidSystem, Region>::
1339equil_OWG(
const Region& reg,
const VSpan& span)
1344 if (! this->oilActive()) {
1345 throw std::invalid_argument {
1346 "Don't know how to interpret EQUIL datum depth in "
1347 "OIL zone in model without active oil phase"
1352 const auto ic =
typename OPress::InitCond {
1353 reg.datum(), reg.pressure()
1356 this->makeOilPressure(ic, reg, span);
1359 if (this->waterActive()) {
1361 const auto ic =
typename WPress::InitCond {
1363 this->oil(reg.zwoc()) - reg.pcowWoc()
1366 this->makeWatPressure(ic, reg, span);
1369 if (this->gasActive()) {
1371 const auto ic =
typename GPress::InitCond {
1373 this->oil(reg.zgoc()) + reg.pcgoGoc()
1375 this->makeGasPressure(ic, reg, span);
1379template <
class Flu
idSystem,
class Region>
1380void PressureTable<FluidSystem, Region>::
1381makeOilPressure(
const typename OPress::InitCond& ic,
1385 const auto drho = OilPressODE {
1386 reg.tempVdTable(), reg.dissolutionCalculator(),
1387 reg.pvtIdx(), this->gravity_
1390 this->oil_ = std::make_unique<OPress>(drho, ic, this->nsample_, span);
1393template <
class Flu
idSystem,
class Region>
1394void PressureTable<FluidSystem, Region>::
1395makeGasPressure(
const typename GPress::InitCond& ic,
1399 const auto drho = GasPressODE {
1400 reg.tempVdTable(), reg.evaporationCalculator(), reg.waterEvaporationCalculator(),
1401 reg.pvtIdx(), this->gravity_
1404 this->gas_ = std::make_unique<GPress>(drho, ic, this->nsample_, span);
1407template <
class Flu
idSystem,
class Region>
1408void PressureTable<FluidSystem, Region>::
1409makeWatPressure(
const typename WPress::InitCond& ic,
1413 const auto drho = WatPressODE {
1414 reg.tempVdTable(), reg.saltVdTable(), reg.pvtIdx(), this->gravity_
1417 this->wat_ = std::make_unique<WPress>(drho, ic, this->nsample_, span);
1422namespace DeckDependent {
1424std::vector<EquilRecord>
1427 const auto& init = state.getInitConfig();
1429 if(!init.hasEquil()) {
1430 throw std::domain_error(
"Deck does not provide equilibration data.");
1433 const auto& equil = init.getEquil();
1434 return { equil.begin(), equil.end() };
1437template<
class Gr
idView>
1440 const GridView& gridview)
1442 std::vector<int> eqlnum(gridview.size(0), 0);
1444 if (eclipseState.fieldProps().has_int(
"EQLNUM")) {
1445 const auto& e = eclipseState.fieldProps().get_int(
"EQLNUM");
1446 std::ranges::transform(e, eqlnum.begin(), [](
int n) { return n - 1; });
1449 const int num_regions = eclipseState.getTableManager().getEqldims().getNumEquilRegions();
1450 if (std::ranges::any_of(eqlnum, [num_regions](
int n){
return n >= num_regions;})) {
1451 throw std::runtime_error(
"Values larger than maximum Equil regions " +
1454 if (std::ranges::any_of(eqlnum, [](
int n){
return n < 0;})) {
1455 throw std::runtime_error(
"zero or negative values provided in EQLNUM");
1462template<
class FluidSystem,
1465 class ElementMapper,
1466 class CartesianIndexMapper>
1467template<
class MaterialLawManager>
1468InitialStateComputer<FluidSystem,
1472 CartesianIndexMapper>::
1473InitialStateComputer(MaterialLawManager& materialLawManager,
1474 const EclipseState& eclipseState,
1476 const GridView& gridView,
1477 const CartesianIndexMapper& cartMapper,
1479 const int num_pressure_points,
1480 const bool applySwatInit)
1481 : temperature_(grid.size(0), eclipseState.getTableManager().rtemp()),
1482 saltConcentration_(grid.size(0)),
1483 saltSaturation_(grid.size(0)),
1484 pp_(FluidSystem::numPhases,
1485 std::vector<Scalar>(grid.size(0))),
1486 sat_(FluidSystem::numPhases,
1487 std::vector<Scalar>(grid.size(0))),
1491 cartesianIndexMapper_(cartMapper),
1492 num_pressure_points_(num_pressure_points)
1495 if (applySwatInit) {
1496 if (eclipseState.fieldProps().has_double(
"SWATINIT")) {
1497 if constexpr (std::is_same_v<Scalar,double>) {
1498 swatInit_ = eclipseState.fieldProps().get_double(
"SWATINIT");
1500 const auto& input = eclipseState.fieldProps().get_double(
"SWATINIT");
1501 swatInit_.resize(input.size());
1502 std::ranges::copy(input, swatInit_.begin());
1509 const auto& num_aquifers = eclipseState.aquifer().numericalAquifers();
1510 updateCellProps_(gridView, num_aquifers);
1513 const std::vector<EquilRecord> rec =
getEquil(eclipseState);
1514 const auto& tables = eclipseState.getTableManager();
1516 const RegionMapping<> eqlmap(
equilnum(eclipseState, grid));
1517 const int invalidRegion = -1;
1518 regionPvtIdx_.resize(rec.size(), invalidRegion);
1519 setRegionPvtIdx(eclipseState, eqlmap);
1522 rsFunc_.reserve(rec.size());
1524 auto getArray = [](
const std::vector<double>& input)
1526 if constexpr (std::is_same_v<Scalar,double>) {
1529 std::vector<Scalar> output;
1530 output.resize(input.size());
1531 std::ranges::copy(input, output.begin());
1536 if (FluidSystem::enableDissolvedGas()) {
1537 for (std::size_t i = 0; i < rec.size(); ++i) {
1538 if (eqlmap.cells(i).empty()) {
1542 const int pvtIdx = regionPvtIdx_[i];
1543 if (!rec[i].liveOilInitConstantRs()) {
1544 const TableContainer& rsvdTables = tables.getRsvdTables();
1545 const TableContainer& pbvdTables = tables.getPbvdTables();
1546 if (rsvdTables.size() > 0) {
1547 const RsvdTable& rsvdTable = rsvdTables.getTable<RsvdTable>(i);
1548 auto depthColumn = getArray(rsvdTable.getColumn(
"DEPTH").vectorCopy());
1549 auto rsColumn = getArray(rsvdTable.getColumn(
"RS").vectorCopy());
1551 depthColumn, rsColumn));
1552 }
else if (pbvdTables.size() > 0) {
1553 const PbvdTable& pbvdTable = pbvdTables.getTable<PbvdTable>(i);
1554 auto depthColumn = getArray(pbvdTable.getColumn(
"DEPTH").vectorCopy());
1555 auto pbubColumn = getArray(pbvdTable.getColumn(
"PBUB").vectorCopy());
1557 depthColumn, pbubColumn));
1560 throw std::runtime_error(
"Cannot initialise: RSVD or PBVD table not available.");
1565 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1566 throw std::runtime_error(
"Cannot initialise: when no explicit RSVD table is given, \n"
1567 "datum depth must be at the gas-oil-contact. "
1568 "In EQUIL region "+
std::to_string(i + 1)+
" (counting from 1), this does not hold.");
1570 const Scalar pContact = rec[i].datumDepthPressure();
1571 const Scalar TContact = 273.15 + 20;
1577 for (std::size_t i = 0; i < rec.size(); ++i) {
1582 rvFunc_.reserve(rec.size());
1583 if (FluidSystem::enableVaporizedOil()) {
1584 for (std::size_t i = 0; i < rec.size(); ++i) {
1585 if (eqlmap.cells(i).empty()) {
1589 const int pvtIdx = regionPvtIdx_[i];
1590 if (!rec[i].wetGasInitConstantRv()) {
1591 const TableContainer& rvvdTables = tables.getRvvdTables();
1592 const TableContainer& pdvdTables = tables.getPdvdTables();
1594 if (rvvdTables.size() > 0) {
1595 const RvvdTable& rvvdTable = rvvdTables.getTable<RvvdTable>(i);
1596 auto depthColumn = getArray(rvvdTable.getColumn(
"DEPTH").vectorCopy());
1597 auto rvColumn = getArray(rvvdTable.getColumn(
"RV").vectorCopy());
1599 depthColumn, rvColumn));
1600 }
else if (pdvdTables.size() > 0) {
1601 const PdvdTable& pdvdTable = pdvdTables.getTable<PdvdTable>(i);
1602 auto depthColumn = getArray(pdvdTable.getColumn(
"DEPTH").vectorCopy());
1603 auto pdewColumn = getArray(pdvdTable.getColumn(
"PDEW").vectorCopy());
1605 depthColumn, pdewColumn));
1607 throw std::runtime_error(
"Cannot initialise: RVVD or PDCD table not available.");
1611 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1612 throw std::runtime_error(
1613 "Cannot initialise: when no explicit RVVD table is given, \n"
1614 "datum depth must be at the gas-oil-contact. "
1615 "In EQUIL region "+
std::to_string(i + 1)+
" (counting from 1), this does not hold.");
1617 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1618 const Scalar TContact = 273.15 + 20;
1624 for (std::size_t i = 0; i < rec.size(); ++i) {
1629 rvwFunc_.reserve(rec.size());
1630 if (FluidSystem::enableVaporizedWater()) {
1631 for (std::size_t i = 0; i < rec.size(); ++i) {
1632 if (eqlmap.cells(i).empty()) {
1636 const int pvtIdx = regionPvtIdx_[i];
1637 if (!rec[i].humidGasInitConstantRvw()) {
1638 const TableContainer& rvwvdTables = tables.getRvwvdTables();
1640 if (rvwvdTables.size() > 0) {
1641 const RvwvdTable& rvwvdTable = rvwvdTables.getTable<RvwvdTable>(i);
1642 auto depthColumn = getArray(rvwvdTable.getColumn(
"DEPTH").vectorCopy());
1643 auto rvwvdColumn = getArray(rvwvdTable.getColumn(
"RVWVD").vectorCopy());
1645 depthColumn, rvwvdColumn));
1647 throw std::runtime_error(
"Cannot initialise: RVWVD table not available.");
1651 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1653 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1655 const auto msg =
"No explicit RVWVD table is given for EQUIL region " +
std::to_string(i + 1) +
". \n"
1656 "and datum depth is not at the gas-oil-contact. \n"
1657 "Rvw is set to 0.0 in all cells. \n";
1658 OpmLog::warning(msg);
1662 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1663 const Scalar TContact = 273.15 + 20;
1670 if (rec[i].waterOilContactDepth() != rec[i].datumDepth()) {
1672 const auto msg =
"No explicit RVWVD table is given for EQUIL region " +
std::to_string(i + 1) +
". \n"
1673 "and datum depth is not at the gas-water-contact. \n"
1674 "Rvw is set to 0.0 in all cells. \n";
1675 OpmLog::warning(msg);
1678 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].waterOilContactCapillaryPressure();
1679 const Scalar TContact = 273.15 + 20;
1687 for (std::size_t i = 0; i < rec.size(); ++i) {
1694 updateInitialTemperature_(eclipseState, eqlmap);
1697 updateInitialSaltConcentration_(eclipseState, eqlmap);
1700 updateInitialSaltSaturation_(eclipseState, eqlmap);
1703 const auto& comm = grid.comm();
1704 calcPressSatRsRv(eqlmap, rec, materialLawManager, gridView, comm, grav);
1707 applyNumericalAquifers_(gridView, num_aquifers, eclipseState.runspec().co2Storage() || eclipseState.runspec().h2Storage());
1713template<
class FluidSystem,
1716 class ElementMapper,
1717 class CartesianIndexMapper>
1723 CartesianIndexMapper>::
1724updateInitialTemperature_(
const EclipseState& eclState,
const RMap& reg)
1726 const int numEquilReg = rsFunc_.size();
1727 tempVdTable_.resize(numEquilReg);
1728 const auto& tables = eclState.getTableManager();
1729 if (!tables.hasTables(
"RTEMPVD")) {
1730 std::vector<Scalar> x = {0.0,1.0};
1731 std::vector<Scalar> y = {
static_cast<Scalar
>(tables.rtemp()),
1732 static_cast<Scalar
>(tables.rtemp())};
1733 for (
auto& table : this->tempVdTable_) {
1734 table.setXYContainers(x, y);
1737 const TableContainer& tempvdTables = tables.getRtempvdTables();
1738 for (std::size_t i = 0; i < tempvdTables.size(); ++i) {
1739 const RtempvdTable& tempvdTable = tempvdTables.getTable<RtempvdTable>(i);
1740 tempVdTable_[i].setXYContainers(tempvdTable.getDepthColumn(), tempvdTable.getTemperatureColumn());
1741 const auto& cells = reg.cells(i);
1742 for (
const auto& cell : cells) {
1743 const Scalar depth = cellCenterDepth_[cell];
1744 this->temperature_[cell] = tempVdTable_[i].eval(depth,
true);
1750template<
class FluidSystem,
1753 class ElementMapper,
1754 class CartesianIndexMapper>
1756void InitialStateComputer<FluidSystem,
1760 CartesianIndexMapper>::
1761updateInitialSaltConcentration_(
const EclipseState& eclState,
const RMap& reg)
1763 const int numEquilReg = rsFunc_.size();
1764 saltVdTable_.resize(numEquilReg);
1765 const auto& tables = eclState.getTableManager();
1766 const TableContainer& saltvdTables = tables.getSaltvdTables();
1769 if (saltvdTables.empty()) {
1770 std::vector<Scalar> x = {0.0,1.0};
1771 std::vector<Scalar> y = {0.0,0.0};
1772 for (
auto& table : this->saltVdTable_) {
1773 table.setXYContainers(x, y);
1776 for (std::size_t i = 0; i < saltvdTables.size(); ++i) {
1777 const SaltvdTable& saltvdTable = saltvdTables.getTable<SaltvdTable>(i);
1778 saltVdTable_[i].setXYContainers(saltvdTable.getDepthColumn(), saltvdTable.getSaltColumn());
1780 const auto& cells = reg.cells(i);
1781 for (
const auto& cell : cells) {
1782 const Scalar depth = cellCenterDepth_[cell];
1783 this->saltConcentration_[cell] = saltVdTable_[i].eval(depth,
true);
1789template<
class FluidSystem,
1792 class ElementMapper,
1793 class CartesianIndexMapper>
1795void InitialStateComputer<FluidSystem,
1799 CartesianIndexMapper>::
1800updateInitialSaltSaturation_(
const EclipseState& eclState,
const RMap& reg)
1802 const int numEquilReg = rsFunc_.size();
1803 saltpVdTable_.resize(numEquilReg);
1804 const auto& tables = eclState.getTableManager();
1805 const TableContainer& saltpvdTables = tables.getSaltpvdTables();
1807 for (std::size_t i = 0; i < saltpvdTables.size(); ++i) {
1808 const SaltpvdTable& saltpvdTable = saltpvdTables.getTable<SaltpvdTable>(i);
1809 saltpVdTable_[i].setXYContainers(saltpvdTable.getDepthColumn(), saltpvdTable.getSaltpColumn());
1811 const auto& cells = reg.cells(i);
1812 for (
const auto& cell : cells) {
1813 const Scalar depth = cellCenterDepth_[cell];
1814 this->saltSaturation_[cell] = saltpVdTable_[i].eval(depth,
true);
1819template<
class FluidSystem,
1822 class ElementMapper,
1823 class CartesianIndexMapper>
1824void InitialStateComputer<FluidSystem,
1828 CartesianIndexMapper>::
1829updateCellProps_(
const GridView& gridView,
1830 const NumericalAquifers& aquifer)
1832 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1833 int numElements = gridView.size(0);
1834 cellCenterDepth_.resize(numElements);
1835 cellCenterXY_.resize(numElements);
1836 cellCorners_.resize(numElements);
1837 cellZSpan_.resize(numElements);
1838 cellZMinMax_.resize(numElements);
1840 auto elemIt = gridView.template begin<0>();
1841 const auto& elemEndIt = gridView.template end<0>();
1842 const auto num_aqu_cells = aquifer.allAquiferCells();
1843 for (; elemIt != elemEndIt; ++elemIt) {
1844 const Element& element = *elemIt;
1845 const unsigned int elemIdx = elemMapper.index(element);
1846 cellCenterDepth_[elemIdx] = Details::cellCenterDepth<Scalar>(element);
1847 cellCenterXY_[elemIdx] = Details::cellCenterXY<Scalar>(element);
1848 cellCorners_[elemIdx] = Details::getCellCornerXY<Scalar>(element);
1849 const auto cartIx = cartesianIndexMapper_.cartesianIndex(elemIdx);
1850 cellZSpan_[elemIdx] = Details::cellZSpan<Scalar>(element);
1851 cellZMinMax_[elemIdx] = Details::cellZMinMax<Scalar>(element);
1852 if (!num_aqu_cells.empty()) {
1853 const auto search = num_aqu_cells.find(cartIx);
1854 if (search != num_aqu_cells.end()) {
1855 const auto* aqu_cell = num_aqu_cells.at(cartIx);
1856 const Scalar depth_change_num_aqu = aqu_cell->depth - cellCenterDepth_[elemIdx];
1857 cellCenterDepth_[elemIdx] += depth_change_num_aqu;
1858 cellZSpan_[elemIdx].first += depth_change_num_aqu;
1859 cellZSpan_[elemIdx].second += depth_change_num_aqu;
1860 cellZMinMax_[elemIdx].first += depth_change_num_aqu;
1861 cellZMinMax_[elemIdx].second += depth_change_num_aqu;
1867template<
class FluidSystem,
1870 class ElementMapper,
1871 class CartesianIndexMapper>
1872void InitialStateComputer<FluidSystem,
1876 CartesianIndexMapper>::
1877applyNumericalAquifers_(
const GridView& gridView,
1878 const NumericalAquifers& aquifer,
1879 const bool co2store_or_h2store)
1881 const auto num_aqu_cells = aquifer.allAquiferCells();
1882 if (num_aqu_cells.empty())
return;
1885 bool oil_as_brine = co2store_or_h2store && FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1886 const auto watPos = oil_as_brine? FluidSystem::oilPhaseIdx : FluidSystem::waterPhaseIdx;
1887 if (!FluidSystem::phaseIsActive(watPos)){
1888 throw std::logic_error {
"Water phase has to be active for numerical aquifer case" };
1891 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1892 auto elemIt = gridView.template begin<0>();
1893 const auto& elemEndIt = gridView.template end<0>();
1894 const auto oilPos = FluidSystem::oilPhaseIdx;
1895 const auto gasPos = FluidSystem::gasPhaseIdx;
1896 for (; elemIt != elemEndIt; ++elemIt) {
1897 const Element& element = *elemIt;
1898 const unsigned int elemIdx = elemMapper.index(element);
1899 const auto cartIx = cartesianIndexMapper_.cartesianIndex(elemIdx);
1900 const auto search = num_aqu_cells.find(cartIx);
1901 if (search != num_aqu_cells.end()) {
1903 this->sat_[watPos][elemIdx] = 1.;
1905 if (!co2store_or_h2store && FluidSystem::phaseIsActive(oilPos)) {
1906 this->sat_[oilPos][elemIdx] = 0.;
1909 if (FluidSystem::phaseIsActive(gasPos)) {
1910 this->sat_[gasPos][elemIdx] = 0.;
1912 const auto* aqu_cell = num_aqu_cells.at(cartIx);
1913 const auto msg = fmt::format(
"FOR AQUIFER CELL AT ({}, {}, {}) OF NUMERICAL "
1914 "AQUIFER {}, WATER SATURATION IS SET TO BE UNITY",
1915 aqu_cell->I+1, aqu_cell->J+1, aqu_cell->K+1, aqu_cell->aquifer_id);
1920 if (aqu_cell->init_pressure) {
1921 const Scalar pres = *(aqu_cell->init_pressure);
1922 this->pp_[watPos][elemIdx] = pres;
1923 if (FluidSystem::phaseIsActive(gasPos)) {
1924 this->pp_[gasPos][elemIdx] = pres;
1926 if (FluidSystem::phaseIsActive(oilPos)) {
1927 this->pp_[oilPos][elemIdx] = pres;
1934template<
class FluidSystem,
1937 class ElementMapper,
1938 class CartesianIndexMapper>
1940void InitialStateComputer<FluidSystem,
1944 CartesianIndexMapper>::
1945setRegionPvtIdx(
const EclipseState& eclState,
const RMap& reg)
1947 const auto& pvtnumData = eclState.fieldProps().get_int(
"PVTNUM");
1949 for (
const auto& r : reg.activeRegions()) {
1950 const auto& cells = reg.cells(r);
1951 regionPvtIdx_[r] = pvtnumData[*cells.begin()] - 1;
1955template<
class FluidSystem,
1958 class ElementMapper,
1959 class CartesianIndexMapper>
1960template<
class RMap,
class MaterialLawManager,
class Comm>
1961void InitialStateComputer<FluidSystem,
1965 CartesianIndexMapper>::
1966calcPressSatRsRv(
const RMap& reg,
1967 const std::vector<EquilRecord>& rec,
1968 MaterialLawManager& materialLawManager,
1969 const GridView& gridView,
1973 using PhaseSat = Details::PhaseSaturations<
1974 MaterialLawManager, FluidSystem, EquilReg<Scalar>,
typename RMap::CellId
1977 auto ptable = Details::PressureTable<FluidSystem, EquilReg<Scalar>>{ grav, this->num_pressure_points_ };
1978 auto psat = PhaseSat { materialLawManager, this->swatInit_ };
1979 auto vspan = std::array<Scalar, 2>{};
1981 std::vector<int> regionIsEmpty(rec.size(), 0);
1982 for (std::size_t r = 0; r < rec.size(); ++r) {
1983 const auto& cells = reg.cells(r);
1987 const auto acc = rec[r].initializationTargetAccuracy();
1991 if (cells.empty()) {
1992 regionIsEmpty[r] = 1;
1995 const auto eqreg = EquilReg {
1996 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
1997 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2000 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2001 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2002 ptable.equilibrate(eqreg, vspan);
2005 this->equilibrateTiltedFaultBlock(cells, eqreg, gridView, acc, ptable, psat);
2007 else if (acc == 0) {
2008 if (cells.empty()) {
2009 regionIsEmpty[r] = 1;
2012 const auto eqreg = EquilReg {
2013 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
2014 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2016 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2017 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2018 ptable.equilibrate(eqreg, vspan);
2020 this->equilibrateCellCentres(cells, eqreg, ptable, psat);
2023 if (cells.empty()) {
2024 regionIsEmpty[r] = 1;
2027 const auto eqreg = EquilReg {
2028 rec[r], this->rsFunc_[r], this->rvFunc_[r], this->rvwFunc_[r],
2029 this->tempVdTable_[r], this->saltVdTable_[r], this->regionPvtIdx_[r]
2031 vspan[0] = std::min(vspan[0], std::min(eqreg.zgoc(), eqreg.zwoc()));
2032 vspan[1] = std::max(vspan[1], std::max(eqreg.zgoc(), eqreg.zwoc()));
2033 ptable.equilibrate(eqreg, vspan);
2035 this->equilibrateHorizontal(cells, eqreg, -acc, ptable, psat);
2038 comm.min(regionIsEmpty.data(),regionIsEmpty.size());
2039 if (comm.rank() == 0) {
2040 for (std::size_t r = 0; r < rec.size(); ++r) {
2041 if (regionIsEmpty[r])
2043 +
" has no active cells");
2048template<
class FluidSystem,
2051 class ElementMapper,
2052 class CartesianIndexMapper>
2053template<
class CellRange,
class EquilibrationMethod>
2054void InitialStateComputer<FluidSystem,
2058 CartesianIndexMapper>::
2059cellLoop(
const CellRange& cells,
2060 EquilibrationMethod&& eqmethod)
2062 const auto oilPos = FluidSystem::oilPhaseIdx;
2063 const auto gasPos = FluidSystem::gasPhaseIdx;
2064 const auto watPos = FluidSystem::waterPhaseIdx;
2066 const auto oilActive = FluidSystem::phaseIsActive(oilPos);
2067 const auto gasActive = FluidSystem::phaseIsActive(gasPos);
2068 const auto watActive = FluidSystem::phaseIsActive(watPos);
2070 auto pressures = Details::PhaseQuantityValue<Scalar>{};
2071 auto saturations = Details::PhaseQuantityValue<Scalar>{};
2076 for (
const auto& cell : cells) {
2077 eqmethod(cell, pressures, saturations, Rs, Rv, Rvw);
2080 this->pp_ [oilPos][cell] = pressures.oil;
2081 this->sat_[oilPos][cell] = saturations.oil;
2085 this->pp_ [gasPos][cell] = pressures.gas;
2086 this->sat_[gasPos][cell] = saturations.gas;
2090 this->pp_ [watPos][cell] = pressures.water;
2091 this->sat_[watPos][cell] = saturations.water;
2094 if (oilActive && gasActive) {
2095 this->rs_[cell] =
Rs;
2096 this->rv_[cell] =
Rv;
2099 if (watActive && gasActive) {
2100 this->rvw_[cell] =
Rvw;
2105template<
class FluidSystem,
2108 class ElementMapper,
2109 class CartesianIndexMapper>
2110template<
class CellRange,
class PressTable,
class PhaseSat>
2111void InitialStateComputer<FluidSystem,
2115 CartesianIndexMapper>::
2116equilibrateCellCentres(
const CellRange& cells,
2117 const EquilReg<Scalar>& eqreg,
2118 const PressTable& ptable,
2121 using CellPos =
typename PhaseSat::Position;
2122 using CellID = std::remove_cv_t<std::remove_reference_t<
2123 decltype(std::declval<CellPos>().cell)>>;
2124 this->cellLoop(cells, [
this, &eqreg, &ptable, &psat]
2126 Details::PhaseQuantityValue<Scalar>& pressures,
2127 Details::PhaseQuantityValue<Scalar>& saturations,
2130 Scalar& Rvw) ->
void
2132 const auto pos = CellPos {
2133 cell, cellCenterDepth_[cell]
2136 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2137 pressures = psat.correctedPhasePressures();
2139 const auto temp = this->temperature_[cell];
2141 Rs = eqreg.dissolutionCalculator()
2142 (pos.depth, pressures.oil, temp, saturations.gas);
2144 Rv = eqreg.evaporationCalculator()
2145 (pos.depth, pressures.gas, temp, saturations.oil);
2147 Rvw = eqreg.waterEvaporationCalculator()
2148 (pos.depth, pressures.gas, temp, saturations.water);
2152template<
class FluidSystem,
2155 class ElementMapper,
2156 class CartesianIndexMapper>
2157template<
class CellRange,
class PressTable,
class PhaseSat>
2158void InitialStateComputer<FluidSystem,
2162 CartesianIndexMapper>::
2163equilibrateHorizontal(
const CellRange& cells,
2164 const EquilReg<Scalar>& eqreg,
2166 const PressTable& ptable,
2169 using CellPos =
typename PhaseSat::Position;
2170 using CellID = std::remove_cv_t<std::remove_reference_t<
2171 decltype(std::declval<CellPos>().cell)>>;
2173 this->cellLoop(cells, [
this, acc, &eqreg, &ptable, &psat]
2175 Details::PhaseQuantityValue<Scalar>& pressures,
2176 Details::PhaseQuantityValue<Scalar>& saturations,
2179 Scalar& Rvw) ->
void
2182 saturations.reset();
2184 Scalar totfrac = 0.0;
2186 const auto pos = CellPos { cell, depth };
2188 saturations.axpy(psat.deriveSaturations(pos, eqreg, ptable), frac);
2189 pressures .axpy(psat.correctedPhasePressures(), frac);
2195 saturations /= totfrac;
2196 pressures /= totfrac;
2199 const auto pos = CellPos {
2200 cell, cellCenterDepth_[cell]
2203 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2204 pressures = psat.correctedPhasePressures();
2207 const auto temp = this->temperature_[cell];
2208 const auto cz = cellCenterDepth_[cell];
2210 Rs = eqreg.dissolutionCalculator()
2211 (cz, pressures.oil, temp, saturations.gas);
2213 Rv = eqreg.evaporationCalculator()
2214 (cz, pressures.gas, temp, saturations.oil);
2216 Rvw = eqreg.waterEvaporationCalculator()
2217 (cz, pressures.gas, temp, saturations.water);
2221template<
class Flu
idSystem,
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper>
2222template<
class CellRange,
class PressTable,
class PhaseSat>
2223void InitialStateComputer<FluidSystem, Grid, GridView, ElementMapper, CartesianIndexMapper>::
2224equilibrateTiltedFaultBlockSimple(
const CellRange& cells,
2225 const EquilReg<Scalar>& eqreg,
2226 const GridView& gridView,
2228 const PressTable& ptable,
2231 using CellPos =
typename PhaseSat::Position;
2232 using CellID = std::remove_cv_t<std::remove_reference_t<
2233 decltype(std::declval<CellPos>().cell)>>;
2235 this->cellLoop(cells, [
this, acc, &eqreg, &ptable, &psat, &gridView]
2237 Details::PhaseQuantityValue<Scalar>& pressures,
2238 Details::PhaseQuantityValue<Scalar>& saturations,
2241 Scalar& Rvw) ->
void
2244 saturations.reset();
2245 Scalar totalWeight = 0.0;
2248 const auto& [zmin, zmax] = cellZMinMax_[cell];
2249 const Scalar cellThickness = zmax - zmin;
2250 const Scalar halfThickness = cellThickness / 2.0;
2253 Scalar dipAngle, dipAzimuth;
2257 std::array<Scalar, 3> referencePoint = {
2258 cellCenterXY_[cell].first,
2259 cellCenterXY_[cell].second,
2260 cellCenterDepth_[cell]
2264 const int numLevelsPerHalf = std::min(20, acc);
2267 std::vector<std::pair<Scalar, Scalar>> levels;
2270 for (
int side = 0; side < 2; ++side) {
2271 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2273 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2275 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2279 Scalar crossSectionWeight = (halfThickness / numLevelsPerHalf);
2282 if (std::abs(dipAngle) > 1e-10) {
2283 crossSectionWeight /= std::cos(dipAngle);
2286 levels.emplace_back(depth, crossSectionWeight);
2290 for (
const auto& [depth, weight] : levels) {
2292 const auto& [x, y] = cellCenterXY_[cell];
2294 depth, x, y, dipAngle, dipAzimuth, referencePoint);
2296 const auto pos = CellPos{cell, tvd};
2298 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2299 auto localPressures = psat.correctedPhasePressures();
2302 saturations.axpy(localSaturations, weight);
2303 pressures.axpy(localPressures, weight);
2304 totalWeight += weight;
2308 if (totalWeight > 1e-10) {
2309 saturations /= totalWeight;
2310 pressures /= totalWeight;
2313 const auto& [x, y] = cellCenterXY_[cell];
2315 cellCenterDepth_[cell], x, y, dipAngle, dipAzimuth, referencePoint);
2316 const auto pos = CellPos{cell, tvdCenter};
2317 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2318 pressures = psat.correctedPhasePressures();
2322 const auto temp = this->temperature_[cell];
2323 const auto& [x, y] = cellCenterXY_[cell];
2325 cellCenterDepth_[cell], x, y, dipAngle, dipAzimuth, referencePoint);
2327 Rs = eqreg.dissolutionCalculator()(tvdCenter, pressures.oil, temp, saturations.gas);
2328 Rv = eqreg.evaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.oil);
2329 Rvw = eqreg.waterEvaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.water);
2333template<
class Flu
idSystem,
class Gr
id,
class Gr
idView,
class ElementMapper,
class CartesianIndexMapper>
2334template<
class CellRange,
class PressTable,
class PhaseSat>
2335void InitialStateComputer<FluidSystem, Grid, GridView, ElementMapper, CartesianIndexMapper>::
2336equilibrateTiltedFaultBlock(
const CellRange& cells,
2337 const EquilReg<Scalar>& eqreg,
2338 const GridView& gridView,
2340 const PressTable& ptable,
2343 using CellPos =
typename PhaseSat::Position;
2344 using CellID = std::remove_cv_t<std::remove_reference_t<
2345 decltype(std::declval<CellPos>().cell)>>;
2347 std::vector<typename GridView::template Codim<0>::Entity> entityMap(gridView.size(0));
2348 for (
const auto& entity : entities(gridView, Dune::Codim<0>())) {
2349 CellID idx = gridView.indexSet().index(entity);
2350 entityMap[idx] = entity;
2354 auto polygonArea = [](
const std::vector<std::array<Scalar, 2>>& pts) {
2355 if (pts.size() < 3)
return Scalar(0);
2357 for (
size_t i = 0; i < pts.size(); ++i) {
2358 size_t j = (i + 1) % pts.size();
2359 area += pts[i][0] * pts[j][1] - pts[j][0] * pts[i][1];
2361 return std::abs(area) * Scalar(0.5);
2365 auto computeCrossSectionArea = [&](
const CellID cell, Scalar depth) -> Scalar {
2367 const auto& entity = entityMap[cell];
2368 const auto& geometry = entity.geometry();
2369 const int numCorners = geometry.corners();
2371 std::vector<std::array<Scalar, 3>> corners(numCorners);
2372 for (
int i = 0; i < numCorners; ++i) {
2373 const auto& corner = geometry.corner(i);
2374 corners[i] = {
static_cast<Scalar
>(corner[0]),
static_cast<Scalar
>(corner[1]),
static_cast<Scalar
>(corner[2])};
2378 std::vector<std::array<Scalar, 2>> intersectionPoints;
2379 const Scalar tol = 1e-10;
2382 for (
size_t i = 0; i < corners.size(); ++i) {
2383 for (
size_t j = i + 1; j < corners.size(); ++j) {
2384 Scalar za = corners[i][2];
2385 Scalar zb = corners[j][2];
2387 if ((za - depth) * (zb - depth) <= 0.0 && std::abs(za - zb) > tol) {
2389 Scalar t = (depth - za) / (zb - za);
2390 Scalar x = corners[i][0] + t * (corners[j][0] - corners[i][0]);
2391 Scalar y = corners[i][1] + t * (corners[j][1] - corners[i][1]);
2392 intersectionPoints.push_back({x, y});
2398 if (intersectionPoints.size() > 1) {
2399 auto pointsEqual = [tol](
const std::array<Scalar, 2>& a,
const std::array<Scalar, 2>& b) {
2400 return std::abs(a[0] - b[0]) < tol && std::abs(a[1] - b[1]) < tol;
2403 intersectionPoints.erase(
2404 std::unique(intersectionPoints.begin(), intersectionPoints.end(), pointsEqual),
2405 intersectionPoints.end()
2409 if (intersectionPoints.size() < 3) {
2415 Scalar cx = 0, cy = 0;
2416 for (
const auto& p : intersectionPoints) {
2417 cx += p[0]; cy += p[1];
2419 cx /= intersectionPoints.size();
2420 cy /= intersectionPoints.size();
2423 auto angleCompare = [cx, cy](
const std::array<Scalar, 2>& a,
const std::array<Scalar, 2>& b) {
2424 return std::atan2(a[1] - cy, a[0] - cx) < std::atan2(b[1] - cy, b[0] - cx);
2427 std::ranges::sort(intersectionPoints, angleCompare);
2429 return polygonArea(intersectionPoints);
2431 }
catch (
const std::exception& e) {
2436 auto cellProcessor = [
this, acc, &eqreg, &ptable, &psat, &computeCrossSectionArea]
2438 Details::PhaseQuantityValue<Scalar>& pressures,
2439 Details::PhaseQuantityValue<Scalar>& saturations,
2442 Scalar&
Rvw) ->
void
2445 saturations.reset();
2446 Scalar totalWeight = 0.0;
2448 const auto& zmin = this->cellZMinMax_[cell].first;
2449 const auto& zmax = this->cellZMinMax_[cell].second;
2450 const Scalar cellThickness = zmax - zmin;
2451 const Scalar halfThickness = cellThickness / 2.0;
2454 Scalar dipAngle, dipAzimuth;
2458 std::array<Scalar, 3> referencePoint = {
2459 this->cellCenterXY_[cell].first,
2460 this->cellCenterXY_[cell].second,
2461 cellCenterDepth_[cell]
2465 const int numLevelsPerHalf = std::min(20, acc);
2468 std::vector<std::pair<Scalar, Scalar>> levels;
2471 for (
int side = 0; side < 2; ++side) {
2472 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2474 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2476 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2479 Scalar crossSectionArea = computeCrossSectionArea(cell, depth);
2482 Scalar weight = crossSectionArea * (halfThickness / numLevelsPerHalf);
2484 levels.emplace_back(depth, weight);
2489 bool hasValidAreas =
false;
2490 for (
const auto& level : levels) {
2491 if (level.second > 1e-10) {
2492 hasValidAreas =
true;
2497 if (!hasValidAreas) {
2500 for (
int side = 0; side < 2; ++side) {
2501 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2502 for (
int i = 0; i < numLevelsPerHalf; ++i) {
2503 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2504 Scalar weight = (halfThickness / numLevelsPerHalf);
2505 if (std::abs(dipAngle) > 1e-10) {
2506 weight /= std::cos(dipAngle);
2508 levels.emplace_back(depth, weight);
2513 for (
const auto& level : levels) {
2514 Scalar depth = level.first;
2515 Scalar weight = level.second;
2518 const auto& xy = this->cellCenterXY_[cell];
2520 depth, xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2522 const auto pos = CellPos{cell, tvd};
2524 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2525 auto localPressures = psat.correctedPhasePressures();
2528 saturations.axpy(localSaturations, weight);
2529 pressures.axpy(localPressures, weight);
2530 totalWeight += weight;
2533 if (totalWeight > 1e-10) {
2534 saturations /= totalWeight;
2535 pressures /= totalWeight;
2538 const auto& xy = this->cellCenterXY_[cell];
2540 this->cellCenterDepth_[cell], xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2541 const auto pos = CellPos{cell, tvdCenter};
2542 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2543 pressures = psat.correctedPhasePressures();
2547 const auto temp = this->temperature_[cell];
2548 const auto& xy = this->cellCenterXY_[cell];
2550 this->cellCenterDepth_[cell], xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2552 Rs = eqreg.dissolutionCalculator()(tvdCenter, pressures.oil, temp, saturations.gas);
2553 Rv = eqreg.evaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.oil);
2554 Rvw = eqreg.waterEvaporationCalculator()(tvdCenter, pressures.gas, temp, saturations.water);
2557 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:472
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:488
Definition: InitStateEquil.hpp:126
Oil(const TabulatedFunction &tempVdTable, const RS &rs, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:424
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:438
Definition: InitStateEquil.hpp:101
Scalar operator()(const Scalar depth, const Scalar press) const
Definition: InitStateEquil_impl.hpp:398
Water(const TabulatedFunction &tempVdTable, const TabulatedFunction &saltVdTable, const int pvtRegionIdx, const Scalar normGrav)
Definition: InitStateEquil_impl.hpp:384
Definition: InitStateEquil.hpp:399
const PhaseQuantityValue< Scalar > & deriveSaturations(const Position &x, const Region ®, const PTable &ptable)
Definition: InitStateEquil_impl.hpp:730
PhaseSaturations(MaterialLawManager &matLawMgr, const std::vector< Scalar > &swatInit)
Definition: InitStateEquil_impl.hpp:706
Definition: InitStateEquil.hpp:180
PressureTable & operator=(const PressureTable &rhs)
Definition: InitStateEquil_impl.hpp:1151
Scalar water(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1231
Scalar gas(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1220
bool waterActive() const
Predicate for whether or not water is an active phase.
Definition: InitStateEquil_impl.hpp:1202
bool gasActive() const
Predicate for whether or not gas is an active phase.
Definition: InitStateEquil_impl.hpp:1195
Scalar oil(const Scalar depth) const
Definition: InitStateEquil_impl.hpp:1210
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:1188
typename FluidSystem::Scalar Scalar
Definition: InitStateEquil.hpp:182
void equilibrate(const Region ®, const VSpan &span)
Definition: InitStateEquil_impl.hpp:1177
PressureTable(const Scalar gravity, const int samplePoints=2000)
Definition: InitStateEquil_impl.hpp:1121
Definition: InitStateEquil.hpp:80
Scalar operator()(const Scalar x) const
Definition: InitStateEquil_impl.hpp:350
RK4IVP(const RHS &f, const std::array< Scalar, 2 > &span, const Scalar y0, const int N)
Definition: InitStateEquil_impl.hpp:315
Definition: EquilibrationHelpers.hpp:135
Definition: EquilibrationHelpers.hpp:216
Definition: EquilibrationHelpers.hpp:269
Definition: EquilibrationHelpers.hpp:162
Definition: EquilibrationHelpers.hpp:322
Definition: EquilibrationHelpers.hpp:376
std::vector< EquilRecord > getEquil(const EclipseState &state)
Definition: InitStateEquil_impl.hpp:1425
std::vector< int > equilnum(const EclipseState &eclipseState, const GridView &gridview)
Definition: InitStateEquil_impl.hpp:1439
std::pair< Scalar, Scalar > cellZMinMax(const Element &element)
Definition: InitStateEquil_impl.hpp:182
Scalar cellCenterDepth(const Element &element)
Definition: InitStateEquil_impl.hpp:129
std::pair< Scalar, Scalar > cellZSpan(const Element &element)
Definition: InitStateEquil_impl.hpp:163
CellCornerData< Scalar > getCellCornerXY(const Element &element)
Definition: InitStateEquil_impl.hpp:252
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:65
std::pair< Scalar, Scalar > cellCenterXY(const Element &element)
Definition: InitStateEquil_impl.hpp:144
Scalar calculateTrueVerticalDepth(Scalar z, Scalar x, Scalar y, Scalar dipAngle, Scalar dipAzimuth, const std::array< Scalar, 3 > &referencePoint)
Definition: InitStateEquil_impl.hpp:276
void subdivisionCentrePoints(const Scalar left, const Scalar right, const int numIntervals, std::vector< std::pair< Scalar, Scalar > > &subdiv)
Definition: InitStateEquil_impl.hpp:90
std::vector< std::pair< Scalar, Scalar > > horizontalSubdivision(const CellID cell, const std::pair< Scalar, Scalar > topbot, const int numIntervals)
Definition: InitStateEquil_impl.hpp:108
void computeBlockDip(const CellCornerData< Scalar > &cellCorners, Scalar &dipAngle, Scalar &dipAzimuth)
Definition: InitStateEquil_impl.hpp:201
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