InitStateEquil_impl.hpp
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 3 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
23#ifndef OPM_INIT_STATE_EQUIL_IMPL_HPP
24#define OPM_INIT_STATE_EQUIL_IMPL_HPP
25
26#include <dune/grid/common/mcmgmapper.hh>
27
28#include <opm/common/OpmLog/OpmLog.hpp>
29
30#include <opm/grid/utility/RegionMapping.hpp>
31
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>
40
41#include <opm/input/eclipse/EclipseState/Tables/SaltpvdTable.hpp>
42
43#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
44#include <opm/material/fluidsystems/BlackOilFluidSystem.hpp>
45
49
50#include <fmt/format.h>
51
52#include <algorithm>
53#include <cassert>
54#include <cstddef>
55#include <limits>
56#include <cmath>
57#include <stdexcept>
58
59namespace Opm {
60namespace EQUIL {
61
62namespace Details {
63
64template <typename CellRange, class Scalar>
65void verticalExtent(const CellRange& cells,
66 const std::vector<std::pair<Scalar, Scalar>>& cellZMinMax,
67 const Parallel::Communication& comm,
68 std::array<Scalar,2>& span)
69{
70 span[0] = std::numeric_limits<Scalar>::max();
71 span[1] = std::numeric_limits<Scalar>::lowest();
72
73 // Define vertical span as
74 //
75 // [minimum(node depth(cells)), maximum(node depth(cells))]
76 //
77 // Note: The implementation of 'RK4IVP<>' implicitly
78 // imposes the requirement that cell centroids are all
79 // within this vertical span. That requirement is not
80 // checked.
81 for (const auto& cell : cells) {
82 if (cellZMinMax[cell].first < span[0]) { span[0] = cellZMinMax[cell].first; }
83 if (cellZMinMax[cell].second > span[1]) { span[1] = cellZMinMax[cell].second; }
84 }
85 span[0] = comm.min(span[0]);
86 span[1] = comm.max(span[1]);
87}
88
89template<class Scalar>
90void subdivisionCentrePoints(const Scalar left,
91 const Scalar right,
92 const int numIntervals,
93 std::vector<std::pair<Scalar, Scalar>>& subdiv)
94{
95 const auto h = (right - left) / numIntervals;
96
97 auto end = left;
98 for (auto i = 0*numIntervals; i < numIntervals; ++i) {
99 const auto start = end;
100 end = left + (i + 1)*h;
101
102 subdiv.emplace_back((start + end) / 2, h);
103 }
104}
105
106template <typename CellID, typename Scalar>
107std::vector<std::pair<Scalar, Scalar>>
108horizontalSubdivision(const CellID cell,
109 const std::pair<Scalar, Scalar> topbot,
110 const int numIntervals)
111{
112 auto subdiv = std::vector<std::pair<Scalar, Scalar>>{};
113 subdiv.reserve(2 * numIntervals);
114
115 if (topbot.first > topbot.second) {
116 throw std::out_of_range {
117 "Negative thickness (inverted top/bottom faces) in cell "
118 + std::to_string(cell)
119 };
120 }
121
122 subdivisionCentrePoints(topbot.first, topbot.second,
123 2*numIntervals, subdiv);
124
125 return subdiv;
126}
127
128template <class Scalar, class Element>
129Scalar cellCenterDepth(const Element& element)
130{
131 typedef typename Element::Geometry Geometry;
132 static constexpr int zCoord = Element::dimension - 1;
133 Scalar zz = 0.0;
134
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];
139
140 return zz/corners;
141}
142
143template <class Scalar, class Element>
144std::pair<Scalar,Scalar> cellCenterXY(const Element& element)
145{
146 typedef typename Element::Geometry Geometry;
147 static constexpr int xCoord = Element::dimension - 3;
148 static constexpr int yCoord = Element::dimension - 2;
149 Scalar yy = 0.0;
150 Scalar xx = 0.0;
151
152
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];
158 }
159 return std::make_pair(xx/corners, yy/corners);
160}
161
162template <class Scalar, class Element>
163std::pair<Scalar,Scalar> cellZSpan(const Element& element)
164{
165 typedef typename Element::Geometry Geometry;
166 static constexpr int zCoord = Element::dimension - 1;
167 Scalar bot = 0.0;
168 Scalar top = 0.0;
169
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];
177
178 return std::make_pair(bot/4, top/4);
179}
180
181template <class Scalar, class Element>
182std::pair<Scalar,Scalar> cellZMinMax(const Element& element)
183{
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();
191
192
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]));
196 }
197 return std::make_pair(min, max);
198}
199
200template<class Scalar>
202 Scalar& dipAngle, Scalar& dipAzimuth)
203{
204 const auto& Xc = cellCorners.X;
205 const auto& Yc = cellCorners.Y;
206 const auto& Zc = cellCorners.Z;
207
208 Scalar v1x = Xc[1] - Xc[0];
209 Scalar v1y = Yc[1] - Yc[0];
210 Scalar v1z = Zc[1] - Zc[0];
211
212 Scalar v2x = Xc[2] - Xc[0];
213 Scalar v2y = Yc[2] - Yc[0];
214 Scalar v2z = Zc[2] - Zc[0];
215
216 // Cross product to get normal vector
217 Scalar nx = v1y * v2z - v1z * v2y;
218 Scalar ny = v1z * v2x - v1x * v2z;
219 Scalar nz = v1x * v2y - v1y * v2x;
220
221 // Normalize the normal vector
222 Scalar norm = std::hypot(nx, ny, nz);
223
224 if (norm > 1e-10) {
225 nx /= norm;
226 ny /= norm;
227 nz /= norm;
228
229 // Dip angle is the angle between normal and vertical (0,0,1)
230 dipAngle = std::acos(std::abs(nz));
231
232 // Dip azimuth (direction of dip)
233 if (std::abs(nx) > 1e-10 || std::abs(ny) > 1e-10) {
234 dipAzimuth = std::atan2(ny, nx);
235 // Convert to 0-2π range
236 dipAzimuth = std::fmod(dipAzimuth + 2*M_PI, 2*M_PI);
237 } else {
238 dipAzimuth = 0.0; // Vertical cell
239 }
240
241 // Clamp dip angle to reasonable values
242 const Scalar maxDip = M_PI/2 - 1e-6;
243 dipAngle = std::min(dipAngle, maxDip);
244 } else {
245 // Degenerate cell - assume horizontal
246 dipAngle = 0.0;
247 dipAzimuth = 0.0;
248 }
249}
250
251template <class Scalar, class Element>
253{
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 {};
264 // Get all 8 corners of the hexahedral cell (maybe expensive)
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];
270 }
271
272 return CellCornerData<Scalar>{X, Y, Z};
273}
274
275template<class Scalar>
276Scalar calculateTrueVerticalDepth(Scalar z, Scalar x, Scalar y,
277 Scalar dipAngle, Scalar dipAzimuth,
278 const std::array<Scalar, 3>& referencePoint)
279{
280 // For True Vertical Depth calculation:
281 // TVD = reference_depth + (z - reference_z) * cos(dipAngle)
282 // + lateral_distance * sin(dipAngle) * cos(azimuth_difference)
283
284 // Calculate lateral displacement from reference point
285 Scalar dx = x - referencePoint[0];
286 Scalar dy = y - referencePoint[1];
287 Scalar dz = z - referencePoint[2];
288
289 // If no dip, TVD is simply the depth
290 if (std::abs(dipAngle) < 1e-10) {
291 return referencePoint[2] + dz;
292 }
293
294 // Calculate the direction from reference point to current point
295 Scalar pointAzimuth = std::atan2(dy, dx);
296
297 // Calculate the angle between dip direction and point direction
298 Scalar azimuthDiff = pointAzimuth - dipAzimuth;
299
300 // Calculate lateral distance
301 Scalar lateralDist = std::hypot(dx, dy);
302
303 // Project lateral distance onto dip direction
304 Scalar lateralInDipDir = lateralDist * std::cos(azimuthDiff);
305
306 // True Vertical Depth calculation
307 // TVD increases with depth (more negative z means deeper)
308 // For a dipping plane: TVD = vertical_component + dip_component
309 Scalar tvd = referencePoint[2] + dz * std::cos(dipAngle) + lateralInDipDir * std::sin(dipAngle);
310
311 return tvd;
312}
313
314template<class Scalar, class RHS>
316 const std::array<Scalar,2>& span,
317 const Scalar y0,
318 const int N)
319 : N_(N)
320 , span_(span)
321{
322 const Scalar h = stepsize();
323 const Scalar h2 = h / 2;
324 const Scalar h6 = h / 6;
325
326 y_.reserve(N + 1);
327 f_.reserve(N + 1);
328
329 y_.push_back(y0);
330 f_.push_back(f(span_[0], y0));
331
332 for (int i = 0; i < N; ++i) {
333 const Scalar x = span_[0] + i*h;
334 const Scalar y = y_.back();
335
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);
340
341 y_.push_back(y + h6*(k1 + 2*(k2 + k3) + k4));
342 f_.push_back(f(x + h, y_.back()));
343 }
344
345 assert (y_.size() == typename std::vector<Scalar>::size_type(N + 1));
346}
347
348template<class Scalar, class RHS>
350operator()(const Scalar x) const
351{
352 // Dense output (O(h**3)) according to Shampine
353 // (Hermite interpolation)
354 const Scalar h = stepsize();
355 int i = (x - span_[0]) / h;
356 const Scalar t = (x - (span_[0] + i*h)) / h;
357
358 // Crude handling of evaluation point outside "span_";
359 if (i < 0) { i = 0; }
360 if (N_ <= i) { i = N_ - 1; }
361
362 const Scalar y0 = y_[i], y1 = y_[i + 1];
363 const Scalar f0 = f_[i], f1 = f_[i + 1];
364
365 Scalar u = (1 - 2*t) * (y1 - y0);
366 u += h * ((t - 1)*f0 + t*f1);
367 u *= t * (t - 1);
368 u += (1 - t)*y0 + t*y1;
369
370 return u;
371}
372
373template<class Scalar, class RHS>
375stepsize() const
376{
377 return (span_[1] - span_[0]) / N_;
378}
379
380namespace PhasePressODE {
381
382template<class FluidSystem>
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)
391 , g_(normGrav)
392{
393}
394
395template<class FluidSystem>
396typename Water<FluidSystem>::Scalar
398operator()(const Scalar depth,
399 const Scalar press) const
400{
401 return this->density(depth, press) * g_;
402}
403
404template<class FluidSystem>
405typename Water<FluidSystem>::Scalar
407density(const Scalar depth,
408 const Scalar press) const
409{
410 // The initializing algorithm can give depths outside the range due to numerical noise i.e. we extrapolate
411 Scalar saltConcentration = saltVdTable_.eval(depth, /*extrapolate=*/true);
412 Scalar temp = tempVdTable_.eval(depth, /*extrapolate=*/true);
413 Scalar rho = FluidSystem::waterPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
414 temp,
415 press,
416 Scalar{0.0} /*=Rsw*/,
417 saltConcentration);
418 rho *= FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
419 return rho;
420}
421
422template<class FluidSystem, class RS>
424Oil(const TabulatedFunction& tempVdTable,
425 const RS& rs,
426 const int pvtRegionIdx,
427 const Scalar normGrav)
428 : tempVdTable_(tempVdTable)
429 , rs_(rs)
430 , pvtRegionIdx_(pvtRegionIdx)
431 , g_(normGrav)
432{
433}
434
435template<class FluidSystem, class RS>
436typename Oil<FluidSystem,RS>::Scalar
438operator()(const Scalar depth,
439 const Scalar press) const
440{
441 return this->density(depth, press) * g_;
442}
443
444template<class FluidSystem, class RS>
445typename Oil<FluidSystem,RS>::Scalar
447density(const Scalar depth,
448 const Scalar press) const
449{
450 const Scalar temp = tempVdTable_.eval(depth, /*extrapolate=*/true);
451 Scalar rs = 0.0;
452 if (FluidSystem::enableDissolvedGas())
453 rs = rs_(depth, press, temp);
454
455 Scalar bOil = 0.0;
456 if (rs >= FluidSystem::oilPvt().saturatedGasDissolutionFactor(pvtRegionIdx_, temp, press)) {
457 bOil = FluidSystem::oilPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
458 }
459 else {
460 bOil = FluidSystem::oilPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rs);
461 }
462 Scalar rho = bOil * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
463 if (FluidSystem::enableDissolvedGas()) {
464 rho += rs * bOil * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
465 }
466
467 return rho;
468}
469
470template<class FluidSystem, class RV, class RVW>
472Gas(const TabulatedFunction& tempVdTable,
473 const RV& rv,
474 const RVW& rvw,
475 const int pvtRegionIdx,
476 const Scalar normGrav)
477 : tempVdTable_(tempVdTable)
478 , rv_(rv)
479 , rvw_(rvw)
480 , pvtRegionIdx_(pvtRegionIdx)
481 , g_(normGrav)
482{
483}
484
485template<class FluidSystem, class RV, class RVW>
486typename Gas<FluidSystem,RV,RVW>::Scalar
488operator()(const Scalar depth,
489 const Scalar press) const
490{
491 return this->density(depth, press) * g_;
492}
493
494template<class FluidSystem, class RV, class RVW>
495typename Gas<FluidSystem,RV,RVW>::Scalar
497density(const Scalar depth,
498 const Scalar press) const
499{
500 const Scalar temp = tempVdTable_.eval(depth, /*extrapolate=*/true);
501 Scalar rv = 0.0;
502 if (FluidSystem::enableVaporizedOil())
503 rv = rv_(depth, press, temp);
504
505 Scalar rvw = 0.0;
506 if (FluidSystem::enableVaporizedWater())
507 rvw = rvw_(depth, press, temp);
508
509 Scalar bGas = 0.0;
510
511 if (FluidSystem::enableVaporizedOil() && FluidSystem::enableVaporizedWater()) {
512 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)
513 && rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press))
514 {
515 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
516 } else {
517 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp, press, rv, rvw);
518 }
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_);
522 return rho;
523 }
524
525 if (FluidSystem::enableVaporizedOil()){
526 if (rv >= FluidSystem::gasPvt().saturatedOilVaporizationFactor(pvtRegionIdx_, temp, press)) {
527 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
528 } else {
529 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
530 temp,
531 press,
532 rv,
533 Scalar{0.0}/*=rvw*/);
534 }
535 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
536 rho += rv * bGas * FluidSystem::referenceDensity(FluidSystem::oilPhaseIdx, pvtRegionIdx_);
537 return rho;
538 }
539
540 if (FluidSystem::enableVaporizedWater()){
541 if (rvw >= FluidSystem::gasPvt().saturatedWaterVaporizationFactor(pvtRegionIdx_, temp, press)) {
542 bGas = FluidSystem::gasPvt().saturatedInverseFormationVolumeFactor(pvtRegionIdx_, temp, press);
543 }
544 else {
545 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_,
546 temp,
547 press,
548 Scalar{0.0} /*=rv*/,
549 rvw);
550 }
551 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
552 rho += rvw * bGas * FluidSystem::referenceDensity(FluidSystem::waterPhaseIdx, pvtRegionIdx_);
553 return rho;
554 }
555
556 // immiscible gas
557 bGas = FluidSystem::gasPvt().inverseFormationVolumeFactor(pvtRegionIdx_, temp,
558 press,
559 Scalar{0.0} /*=rv*/,
560 Scalar{0.0} /*=rvw*/);
561 Scalar rho = bGas * FluidSystem::referenceDensity(FluidSystem::gasPhaseIdx, pvtRegionIdx_);
562
563 return rho;
564}
565
566}
567
568template<class FluidSystem, class Region>
569template<class ODE>
570PressureTable<FluidSystem,Region>::
571PressureFunction<ODE>::PressureFunction(const ODE& ode,
572 const InitCond& ic,
573 const int nsample,
574 const VSpan& span)
575 : initial_(ic)
576{
577 this->value_[Direction::Up] = std::make_unique<Distribution>
578 (ode, VSpan {{ ic.depth, span[0] }}, ic.pressure, nsample);
579
580 this->value_[Direction::Down] = std::make_unique<Distribution>
581 (ode, VSpan {{ ic.depth, span[1] }}, ic.pressure, nsample);
582}
583
584template<class FluidSystem, class Region>
585template<class ODE>
586PressureTable<FluidSystem,Region>::
587PressureFunction<ODE>::PressureFunction(const PressureFunction& rhs)
588 : initial_(rhs.initial_)
589{
590 this->value_[Direction::Up] =
591 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
592
593 this->value_[Direction::Down] =
594 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
595}
596
597template<class FluidSystem, class Region>
598template<class ODE>
599typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
602operator=(const PressureFunction& rhs)
603{
604 this->initial_ = rhs.initial_;
605
606 this->value_[Direction::Up] =
607 std::make_unique<Distribution>(*rhs.value_[Direction::Up]);
608
609 this->value_[Direction::Down] =
610 std::make_unique<Distribution>(*rhs.value_[Direction::Down]);
611
612 return *this;
613}
614
615template<class FluidSystem, class Region>
616template<class ODE>
617typename PressureTable<FluidSystem,Region>::template PressureFunction<ODE>&
620operator=(PressureFunction&& rhs)
621{
622 this->initial_ = rhs.initial_;
623 this->value_ = std::move(rhs.value_);
624
625 return *this;
626}
627
628template<class FluidSystem, class Region>
629template<class ODE>
631PressureTable<FluidSystem,Region>::
632PressureFunction<ODE>::
633value(const Scalar depth) const
634{
635 if (depth < this->initial_.depth) {
636 // Value above initial condition depth.
637 return (*this->value_[Direction::Up])(depth);
638 }
639 else if (depth > this->initial_.depth) {
640 // Value below initial condition depth.
641 return (*this->value_[Direction::Down])(depth);
642 }
643 else {
644 // Value *at* initial condition depth.
645 return this->initial_.pressure;
646 }
647}
648
649
650template<class FluidSystem, class Region>
651template<typename PressFunc>
652void PressureTable<FluidSystem,Region>::
653checkPtr(const PressFunc* phasePress,
654 const std::string& phaseName) const
655{
656 if (phasePress != nullptr) { return; }
657
658 throw std::invalid_argument {
659 "Phase pressure function for \"" + phaseName
660 + "\" most not be null"
661 };
662}
663
664template<class FluidSystem, class Region>
665typename PressureTable<FluidSystem,Region>::Strategy
666PressureTable<FluidSystem,Region>::
667selectEquilibrationStrategy(const Region& reg) const
668{
669 if (!this->oilActive()) {
670 if (reg.datum() > reg.zwoc()) { // Datum in water zone
671 return &PressureTable::equil_WOG;
672 }
673 return &PressureTable::equil_GOW;
674 }
675
676 if (reg.datum() > reg.zwoc()) { // Datum in water zone
677 return &PressureTable::equil_WOG;
678 }
679 else if (reg.datum() < reg.zgoc()) { // Datum in gas zone
680 return &PressureTable::equil_GOW;
681 }
682 else { // Datum in oil zone
683 return &PressureTable::equil_OWG;
684 }
685}
686
687template<class FluidSystem, class Region>
688void PressureTable<FluidSystem,Region>::
689copyInPointers(const PressureTable& rhs)
690{
691 if (rhs.oil_ != nullptr) {
692 this->oil_ = std::make_unique<OPress>(*rhs.oil_);
693 }
694
695 if (rhs.gas_ != nullptr) {
696 this->gas_ = std::make_unique<GPress>(*rhs.gas_);
697 }
698
699 if (rhs.wat_ != nullptr) {
700 this->wat_ = std::make_unique<WPress>(*rhs.wat_);
701 }
702}
703
704template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
706PhaseSaturations(MaterialLawManager& matLawMgr,
707 const std::vector<Scalar>& swatInit)
708 : matLawMgr_(matLawMgr)
709 , swatInit_ (swatInit)
710{
711}
712
713template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
716 : matLawMgr_(rhs.matLawMgr_)
717 , swatInit_ (rhs.swatInit_)
718 , sat_ (rhs.sat_)
719 , press_ (rhs.press_)
720{
721 // Note: We don't need to do anything to the 'fluidState_' here.
722 this->setEvaluationPoint(*rhs.evalPt_.position,
723 *rhs.evalPt_.region,
724 *rhs.evalPt_.ptable);
725}
726
727template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
731 const Region& reg,
732 const PTable& ptable)
733{
734 this->setEvaluationPoint(x, reg, ptable);
735 this->initializePhaseQuantities();
736
737 if (ptable.gasActive()) { this->deriveGasSat(); }
738
739 if (ptable.waterActive()) { this->deriveWaterSat(); }
740
741
742 if (this->isOverlappingTransition()) {
743 this->fixUnphysicalTransition();
744 }
745
746 if (ptable.oilActive()) { this->deriveOilSat(); }
747
748 this->accountForScaledSaturations();
749
750 return this->sat_;
751}
752
753template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
755setEvaluationPoint(const Position& x,
756 const Region& reg,
757 const PTable& ptable)
758{
759 this->evalPt_.position = &x;
760 this->evalPt_.region = &reg;
761 this->evalPt_.ptable = &ptable;
762}
763
764template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
765void PhaseSaturations<MaterialLawManager,FluidSystem,Region,CellID>::
766initializePhaseQuantities()
767{
768 this->sat_.reset();
769 this->press_.reset();
770
771 const auto depth = this->evalPt_.position->depth;
772 const auto& ptable = *this->evalPt_.ptable;
773
774 if (ptable.oilActive()) {
775 this->press_.oil = ptable.oil(depth);
776 }
777
778 if (ptable.gasActive()) {
779 this->press_.gas = ptable.gas(depth);
780 }
781
782 if (ptable.waterActive()) {
783 this->press_.water = ptable.water(depth);
784 }
785}
786
787template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
788void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveOilSat()
789{
790 this->sat_.oil = 1.0 - this->sat_.water - this->sat_.gas;
791}
792
793template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
794void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveGasSat()
795{
796 auto& sg = this->sat_.gas;
797
798 const auto isIncr = true; // dPcgo/dSg >= 0 for all Sg.
799 const auto oilActive = this->evalPt_.ptable->oilActive();
800
801 if (this->isConstCapPress(this->gasPos())) {
802 // Sharp interface between phases. Can derive phase saturation
803 // directly from knowing where 'depth' of evaluation point is
804 // relative to depth of O/G contact.
805 const auto gas_contact = oilActive? this->evalPt_.region->zgoc() : this->evalPt_.region->zwoc();
806 sg = this->fromDepthTable(gas_contact,
807 this->gasPos(), isIncr);
808 }
809 else {
810 // Capillary pressure curve is non-constant, meaning there is a
811 // transition zone between the gas and oil phases. Invert capillary
812 // pressure relation
813 //
814 // Pcgo(Sg) = Pg - Po
815 //
816 // Note that Pcgo is defined to be (Pg - Po), not (Po - Pg).
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);
820 }
821}
822
823template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
824void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::deriveWaterSat()
825{
826 auto& sw = this->sat_.water;
827
828 const auto oilActive = this->evalPt_.ptable->oilActive();
829 if (!oilActive) {
830 // for 2p gas+water we set the water saturation to 1.0 - sg
831 sw = 1.0 - this->sat_.gas;
832 }
833 else {
834 const auto isIncr = false; // dPcow/dSw <= 0 for all Sw.
835
836 if (this->isConstCapPress(this->waterPos())) {
837 // Sharp interface between phases. Can derive phase saturation
838 // directly from knowing where 'depth' of evaluation point is
839 // relative to depth of O/W contact.
840 sw = this->fromDepthTable(this->evalPt_.region->zwoc(),
841 this->waterPos(), isIncr);
842 }
843 else {
844 // Capillary pressure curve is non-constant, meaning there is a
845 // transition zone between the oil and water phases. Invert
846 // capillary pressure relation
847 //
848 // Pcow(Sw) = Po - Pw
849 //
850 // unless the model uses "SWATINIT". In the latter case, pick the
851 // saturation directly from the SWATINIT array of the pertinent
852 // cell.
853 const auto pcow = this->press_.oil - this->press_.water;
854
855 if (this->swatInit_.empty()) {
856 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
857 }
858 else {
859 auto [swout, newSwatInit] = this->applySwatInit(pcow);
860 if (newSwatInit)
861 sw = this->invertCapPress(pcow, this->waterPos(), isIncr);
862 else {
863 sw = swout;
864 }
865 }
866 }
867 }
868}
869
870template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
871void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
872fixUnphysicalTransition()
873{
874 auto& sg = this->sat_.gas;
875 auto& sw = this->sat_.water;
876
877 // Overlapping gas/oil and oil/water transition zones can lead to
878 // unphysical phase saturations when individual saturations are derived
879 // directly from inverting O/G and O/W capillary pressure curves.
880 //
881 // Recalculate phase saturations using the implied gas/water capillary
882 // pressure: Pg - Pw.
883 const auto pcgw = this->press_.gas - this->press_.water;
884 if (! this->swatInit_.empty()) {
885 // Re-scale Pc to reflect imposed sw for vanishing oil phase. This
886 // seems consistent with ECLIPSE, but fails to honour SWATINIT in
887 // case of non-trivial gas/oil capillary pressure.
888 auto [swout, newSwatInit] = this->applySwatInit(pcgw, sw);
889 if (newSwatInit){
890 const auto isIncr = false; // dPcow/dSw <= 0 for all Sw.
891 sw = this->invertCapPress(pcgw, this->waterPos(), isIncr);
892 }
893 else {
894 sw = swout;
895 }
896 }
897
898 sw = satFromSumOfPcs<FluidSystem>
899 (this->matLawMgr_, this->waterPos(), this->gasPos(),
900 this->evalPt_.position->cell, pcgw);
901 sg = 1.0 - sw;
902
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);
907
908 // Pcgo = Pg - Po => Po = Pg - Pcgo
909 this->computeMaterialLawCapPress();
910 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
911}
912
913template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
914void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
915accountForScaledSaturations()
916{
917 const auto gasActive = this->evalPt_.ptable->gasActive();
918 const auto watActive = this->evalPt_.ptable->waterActive();
919 const auto oilActive = this->evalPt_.ptable->oilActive();
920
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;
924
925 this->fluidState_.setSaturation(this->waterPos(), sw);
926 this->fluidState_.setSaturation(this->oilPos(), so);
927 this->fluidState_.setSaturation(this->gasPos(), sg);
928
929 const auto& scaledDrainageInfo = this->matLawMgr_
930 .oilWaterScaledEpsInfoDrainage(this->evalPt_.position->cell);
931
932 const auto thresholdSat = 1.0e-6;
933 if (watActive && ((sw + thresholdSat) > scaledDrainageInfo.Swu)) {
934 // Water saturation exceeds maximum possible value. Reset oil phase
935 // pressure to that which corresponds to maximum possible water
936 // saturation value.
937 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swu);
938 if (oilActive) {
939 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swu);
940 } else if (gasActive) {
941 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swu);
942 }
943 sw = scaledDrainageInfo.Swu;
944 this->computeMaterialLawCapPress();
945
946 if (oilActive) {
947 // Pcow = Po - Pw => Po = Pw + Pcow
948 this->press_.oil = this->press_.water + this->materialLawCapPressOilWater();
949 } else {
950 // Pcgw = Pg - Pw => Pg = Pw + Pcgw
951 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
952 }
953
954 }
955 if (gasActive && ((sg + thresholdSat) > scaledDrainageInfo.Sgu)) {
956 // Gas saturation exceeds maximum possible value. Reset oil phase
957 // pressure to that which corresponds to maximum possible gas
958 // saturation value.
959 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgu);
960 if (oilActive) {
961 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgu);
962 } else if (watActive) {
963 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgu);
964 }
965 sg = scaledDrainageInfo.Sgu;
966 this->computeMaterialLawCapPress();
967
968 if (oilActive) {
969 // Pcgo = Pg - Po => Po = Pg - Pcgo
970 this->press_.oil = this->press_.gas - this->materialLawCapPressGasOil();
971 } else {
972 // Pcgw = Pg - Pw => Pw = Pg - Pcgw
973 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
974 }
975 }
976
977 if (watActive && ((sw - thresholdSat) < scaledDrainageInfo.Swl)) {
978 // Water saturation less than minimum possible value in cell. Reset
979 // water phase pressure to that which corresponds to minimum
980 // possible water saturation value.
981 this->fluidState_.setSaturation(this->waterPos(), scaledDrainageInfo.Swl);
982 if (oilActive) {
983 this->fluidState_.setSaturation(this->oilPos(), so + sw - scaledDrainageInfo.Swl);
984 } else if (gasActive) {
985 this->fluidState_.setSaturation(this->gasPos(), sg + sw - scaledDrainageInfo.Swl);
986 }
987 sw = scaledDrainageInfo.Swl;
988 this->computeMaterialLawCapPress();
989
990 if (oilActive) {
991 // Pcwo = Po - Pw => Pw = Po - Pcow
992 this->press_.water = this->press_.oil - this->materialLawCapPressOilWater();
993 } else {
994 // Pcgw = Pg - Pw => Pw = Pg - Pcgw
995 this->press_.water = this->press_.gas - this->materialLawCapPressGasWater();
996 }
997 }
998
999 if (gasActive && ((sg - thresholdSat) < scaledDrainageInfo.Sgl)) {
1000 // Gas saturation less than minimum possible value in cell. Reset
1001 // gas phase pressure to that which corresponds to minimum possible
1002 // gas saturation.
1003 this->fluidState_.setSaturation(this->gasPos(), scaledDrainageInfo.Sgl);
1004 if (oilActive) {
1005 this->fluidState_.setSaturation(this->oilPos(), so + sg - scaledDrainageInfo.Sgl);
1006 } else if (watActive) {
1007 this->fluidState_.setSaturation(this->waterPos(), sw + sg - scaledDrainageInfo.Sgl);
1008 }
1009 sg = scaledDrainageInfo.Sgl;
1010 this->computeMaterialLawCapPress();
1011
1012 if (oilActive) {
1013 // Pcgo = Pg - Po => Pg = Po + Pcgo
1014 this->press_.gas = this->press_.oil + this->materialLawCapPressGasOil();
1015 } else {
1016 // Pcgw = Pg - Pw => Pg = Pw + Pcgw
1017 this->press_.gas = this->press_.water + this->materialLawCapPressGasWater();
1018 }
1019 }
1020}
1021
1022template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1023std::pair<typename FluidSystem::Scalar, bool>
1024PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1025applySwatInit(const Scalar pcow)
1026{
1027 return this->applySwatInit(pcow, this->swatInit_[this->evalPt_.position->cell]);
1028}
1029
1030template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1031std::pair<typename FluidSystem::Scalar, bool>
1032PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1033applySwatInit(const Scalar pcow, const Scalar sw)
1034{
1035 return this->matLawMgr_.applySwatinit(this->evalPt_.position->cell, pcow, sw);
1036}
1037
1038template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1039void PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1040computeMaterialLawCapPress()
1041{
1042 const auto& matParams = this->matLawMgr_
1043 .materialLawParams(this->evalPt_.position->cell);
1044
1045 this->matLawCapPress_.fill(0.0);
1046 MaterialLaw::capillaryPressures(this->matLawCapPress_,
1047 matParams, this->fluidState_);
1048}
1049
1050template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1051typename FluidSystem::Scalar
1052PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1053materialLawCapPressGasOil() const
1054{
1055 return this->matLawCapPress_[this->oilPos()]
1056 + this->matLawCapPress_[this->gasPos()];
1057}
1058
1059template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1060typename FluidSystem::Scalar
1061PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1062materialLawCapPressOilWater() const
1063{
1064 return this->matLawCapPress_[this->oilPos()]
1065 - this->matLawCapPress_[this->waterPos()];
1066}
1067
1068template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1069typename FluidSystem::Scalar
1070PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1071materialLawCapPressGasWater() const
1072{
1073 return this->matLawCapPress_[this->gasPos()]
1074 - this->matLawCapPress_[this->waterPos()];
1075}
1076
1077template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1078bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1079isConstCapPress(const PhaseIdx phaseIdx) const
1080{
1081 return isConstPc<FluidSystem>
1082 (this->matLawMgr_, phaseIdx, this->evalPt_.position->cell);
1083}
1084
1085template <class MaterialLawManager, class FluidSystem, class Region, typename CellID>
1086bool PhaseSaturations<MaterialLawManager, FluidSystem, Region, CellID>::
1087isOverlappingTransition() const
1088{
1089 return this->evalPt_.ptable->gasActive()
1090 && this->evalPt_.ptable->waterActive()
1091 && ((this->sat_.gas + this->sat_.water) > 1.0);
1092}
1093
1094template <class MaterialLawManager, class FluidSystem, 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
1100{
1101 return satFromDepth<FluidSystem>
1102 (this->matLawMgr_, this->evalPt_.position->depth,
1103 contactdepth, static_cast<int>(phasePos),
1104 this->evalPt_.position->cell, isincr);
1105}
1106
1107template <class MaterialLawManager, class FluidSystem, 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
1113{
1114 return satFromPc<FluidSystem>
1115 (this->matLawMgr_, static_cast<int>(phasePos),
1116 this->evalPt_.position->cell, pc, isincr);
1117}
1118
1119template<class FluidSystem, class Region>
1121PressureTable(const Scalar gravity,
1122 const int samplePoints)
1123 : gravity_(gravity)
1124 , nsample_(samplePoints)
1125{
1126}
1127
1128template <class FluidSystem, class Region>
1131 : gravity_(rhs.gravity_)
1132 , nsample_(rhs.nsample_)
1133{
1134 this->copyInPointers(rhs);
1135}
1136
1137template <class FluidSystem, 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_))
1145{
1146}
1147
1148template <class FluidSystem, class Region>
1152{
1153 this->gravity_ = rhs.gravity_;
1154 this->nsample_ = rhs.nsample_;
1155 this->copyInPointers(rhs);
1156
1157 return *this;
1158}
1159
1160template <class FluidSystem, class Region>
1164{
1165 this->gravity_ = rhs.gravity_;
1166 this->nsample_ = rhs.nsample_;
1167
1168 this->oil_ = std::move(rhs.oil_);
1169 this->gas_ = std::move(rhs.gas_);
1170 this->wat_ = std::move(rhs.wat_);
1171
1172 return *this;
1173}
1174
1175template <class FluidSystem, class Region>
1177equilibrate(const Region& reg,
1178 const VSpan& span)
1179{
1180 // One of the PressureTable::equil_*() member functions.
1181 auto equil = this->selectEquilibrationStrategy(reg);
1182
1183 (this->*equil)(reg, span);
1184}
1185
1186template <class FluidSystem, class Region>
1188oilActive() const
1189{
1190 return FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1191}
1192
1193template <class FluidSystem, class Region>
1195gasActive() const
1196{
1197 return FluidSystem::phaseIsActive(FluidSystem::gasPhaseIdx);
1198}
1199
1200template <class FluidSystem, class Region>
1202waterActive() const
1203{
1204 return FluidSystem::phaseIsActive(FluidSystem::waterPhaseIdx);
1205}
1206
1207template <class FluidSystem, class Region>
1208typename FluidSystem::Scalar
1210oil(const Scalar depth) const
1211{
1212 this->checkPtr(this->oil_.get(), "OIL");
1213
1214 return this->oil_->value(depth);
1215}
1216
1217template <class FluidSystem, class Region>
1218typename FluidSystem::Scalar
1220gas(const Scalar depth) const
1221{
1222 this->checkPtr(this->gas_.get(), "GAS");
1223
1224 return this->gas_->value(depth);
1225}
1226
1227
1228template <class FluidSystem, class Region>
1229typename FluidSystem::Scalar
1231water(const Scalar depth) const
1232{
1233 this->checkPtr(this->wat_.get(), "WATER");
1234
1235 return this->wat_->value(depth);
1236}
1237
1238template <class FluidSystem, class Region>
1240equil_WOG(const Region& reg, const VSpan& span)
1241{
1242 // Datum depth in water zone. Calculate phase pressure for water first,
1243 // followed by oil and gas if applicable.
1244
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"
1249 };
1250 }
1251
1252 {
1253 const auto ic = typename WPress::InitCond {
1254 reg.datum(), reg.pressure()
1255 };
1256
1257 this->makeWatPressure(ic, reg, span);
1258 }
1259
1260 if (this->oilActive()) {
1261 // Pcow = Po - Pw => Po = Pw + Pcow
1262 const auto ic = typename OPress::InitCond {
1263 reg.zwoc(),
1264 this->water(reg.zwoc()) + reg.pcowWoc()
1265 };
1266
1267 this->makeOilPressure(ic, reg, span);
1268 }
1269
1270 if (this->gasActive() && this->oilActive()) {
1271 // Pcgo = Pg - Po => Pg = Po + Pcgo
1272 const auto ic = typename GPress::InitCond {
1273 reg.zgoc(),
1274 this->oil(reg.zgoc()) + reg.pcgoGoc()
1275 };
1276
1277 this->makeGasPressure(ic, reg, span);
1278 } else if (this->gasActive() && !this->oilActive()) {
1279 // No oil phase set Pg = Pw + Pcgw
1280 const auto ic = typename GPress::InitCond {
1281 reg.zwoc(), // The WOC is really the GWC for gas/water cases
1282 this->water(reg.zwoc()) + reg.pcowWoc() // Pcow(WOC) is really Pcgw(GWC) for gas/water cases
1283 };
1284 this->makeGasPressure(ic, reg, span);
1285 }
1286}
1287
1288template <class FluidSystem, class Region>
1289void PressureTable<FluidSystem, Region>::
1290equil_GOW(const Region& reg, const VSpan& span)
1291{
1292 // Datum depth in gas zone. Calculate phase pressure for gas first,
1293 // followed by oil and water if applicable.
1294
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"
1299 };
1300 }
1301
1302 {
1303 const auto ic = typename GPress::InitCond {
1304 reg.datum(), reg.pressure()
1305 };
1306
1307 this->makeGasPressure(ic, reg, span);
1308 }
1309
1310 if (this->oilActive()) {
1311 // Pcgo = Pg - Po => Po = Pg - Pcgo
1312 const auto ic = typename OPress::InitCond {
1313 reg.zgoc(),
1314 this->gas(reg.zgoc()) - reg.pcgoGoc()
1315 };
1316 this->makeOilPressure(ic, reg, span);
1317 }
1318
1319 if (this->waterActive() && this->oilActive()) {
1320 // Pcow = Po - Pw => Pw = Po - Pcow
1321 const auto ic = typename WPress::InitCond {
1322 reg.zwoc(),
1323 this->oil(reg.zwoc()) - reg.pcowWoc()
1324 };
1325
1326 this->makeWatPressure(ic, reg, span);
1327 } else if (this->waterActive() && !this->oilActive()) {
1328 // No oil phase set Pw = Pg - Pcgw
1329 const auto ic = typename WPress::InitCond {
1330 reg.zwoc(), // The WOC is really the GWC for gas/water cases
1331 this->gas(reg.zwoc()) - reg.pcowWoc() // Pcow(WOC) is really Pcgw(GWC) for gas/water cases
1332 };
1333 this->makeWatPressure(ic, reg, span);
1334 }
1335}
1336
1337template <class FluidSystem, class Region>
1338void PressureTable<FluidSystem, Region>::
1339equil_OWG(const Region& reg, const VSpan& span)
1340{
1341 // Datum depth in oil zone. Calculate phase pressure for oil first,
1342 // followed by gas and water if applicable.
1343
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"
1348 };
1349 }
1350
1351 {
1352 const auto ic = typename OPress::InitCond {
1353 reg.datum(), reg.pressure()
1354 };
1355
1356 this->makeOilPressure(ic, reg, span);
1357 }
1358
1359 if (this->waterActive()) {
1360 // Pcow = Po - Pw => Pw = Po - Pcow
1361 const auto ic = typename WPress::InitCond {
1362 reg.zwoc(),
1363 this->oil(reg.zwoc()) - reg.pcowWoc()
1364 };
1365
1366 this->makeWatPressure(ic, reg, span);
1367 }
1368
1369 if (this->gasActive()) {
1370 // Pcgo = Pg - Po => Pg = Po + Pcgo
1371 const auto ic = typename GPress::InitCond {
1372 reg.zgoc(),
1373 this->oil(reg.zgoc()) + reg.pcgoGoc()
1374 };
1375 this->makeGasPressure(ic, reg, span);
1376 }
1377}
1378
1379template <class FluidSystem, class Region>
1380void PressureTable<FluidSystem, Region>::
1381makeOilPressure(const typename OPress::InitCond& ic,
1382 const Region& reg,
1383 const VSpan& span)
1384{
1385 const auto drho = OilPressODE {
1386 reg.tempVdTable(), reg.dissolutionCalculator(),
1387 reg.pvtIdx(), this->gravity_
1388 };
1389
1390 this->oil_ = std::make_unique<OPress>(drho, ic, this->nsample_, span);
1391}
1392
1393template <class FluidSystem, class Region>
1394void PressureTable<FluidSystem, Region>::
1395makeGasPressure(const typename GPress::InitCond& ic,
1396 const Region& reg,
1397 const VSpan& span)
1398{
1399 const auto drho = GasPressODE {
1400 reg.tempVdTable(), reg.evaporationCalculator(), reg.waterEvaporationCalculator(),
1401 reg.pvtIdx(), this->gravity_
1402 };
1403
1404 this->gas_ = std::make_unique<GPress>(drho, ic, this->nsample_, span);
1405}
1406
1407template <class FluidSystem, class Region>
1408void PressureTable<FluidSystem, Region>::
1409makeWatPressure(const typename WPress::InitCond& ic,
1410 const Region& reg,
1411 const VSpan& span)
1412{
1413 const auto drho = WatPressODE {
1414 reg.tempVdTable(), reg.saltVdTable(), reg.pvtIdx(), this->gravity_
1415 };
1416
1417 this->wat_ = std::make_unique<WPress>(drho, ic, this->nsample_, span);
1418}
1419
1420}
1421
1422namespace DeckDependent {
1423
1424std::vector<EquilRecord>
1425getEquil(const EclipseState& state)
1426{
1427 const auto& init = state.getInitConfig();
1428
1429 if(!init.hasEquil()) {
1430 throw std::domain_error("Deck does not provide equilibration data.");
1431 }
1432
1433 const auto& equil = init.getEquil();
1434 return { equil.begin(), equil.end() };
1435}
1436
1437template<class GridView>
1438std::vector<int>
1439equilnum(const EclipseState& eclipseState,
1440 const GridView& gridview)
1441{
1442 std::vector<int> eqlnum(gridview.size(0), 0);
1443
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; });
1447 }
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 " +
1452 std::to_string(num_regions) + " provided in EQLNUM");
1453 }
1454 if (std::ranges::any_of(eqlnum, [](int n){return n < 0;})) {
1455 throw std::runtime_error("zero or negative values provided in EQLNUM");
1456 }
1457 OPM_END_PARALLEL_TRY_CATCH("Invalied EQLNUM numbers: ", gridview.comm());
1458
1459 return eqlnum;
1460}
1461
1462template<class FluidSystem,
1463 class Grid,
1464 class GridView,
1465 class ElementMapper,
1466 class CartesianIndexMapper>
1467template<class MaterialLawManager>
1468InitialStateComputer<FluidSystem,
1469 Grid,
1470 GridView,
1471 ElementMapper,
1472 CartesianIndexMapper>::
1473InitialStateComputer(MaterialLawManager& materialLawManager,
1474 const EclipseState& eclipseState,
1475 const Grid& grid,
1476 const GridView& gridView,
1477 const CartesianIndexMapper& cartMapper,
1478 const Scalar grav,
1479 const int num_pressure_points,
1480 const bool applySwatInit)
1481 : temperature_(grid.size(/*codim=*/0), eclipseState.getTableManager().rtemp()),
1482 saltConcentration_(grid.size(/*codim=*/0)),
1483 saltSaturation_(grid.size(/*codim=*/0)),
1484 pp_(FluidSystem::numPhases,
1485 std::vector<Scalar>(grid.size(/*codim=*/0))),
1486 sat_(FluidSystem::numPhases,
1487 std::vector<Scalar>(grid.size(/*codim=*/0))),
1488 rs_(grid.size(/*codim=*/0)),
1489 rv_(grid.size(/*codim=*/0)),
1490 rvw_(grid.size(/*codim=*/0)),
1491 cartesianIndexMapper_(cartMapper),
1492 num_pressure_points_(num_pressure_points)
1493{
1494 //Check for presence of kw SWATINIT
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");
1499 } else {
1500 const auto& input = eclipseState.fieldProps().get_double("SWATINIT");
1501 swatInit_.resize(input.size());
1502 std::ranges::copy(input, swatInit_.begin());
1503 }
1504 }
1505 }
1506
1507 // Querry cell depth, cell top-bottom.
1508 // numerical aquifer cells might be specified with different depths.
1509 const auto& num_aquifers = eclipseState.aquifer().numericalAquifers();
1510 updateCellProps_(gridView, num_aquifers);
1511
1512 // Get the equilibration records.
1513 const std::vector<EquilRecord> rec = getEquil(eclipseState);
1514 const auto& tables = eclipseState.getTableManager();
1515 // Create (inverse) region mapping.
1516 const RegionMapping<> eqlmap(equilnum(eclipseState, grid));
1517 const int invalidRegion = -1;
1518 regionPvtIdx_.resize(rec.size(), invalidRegion);
1519 setRegionPvtIdx(eclipseState, eqlmap);
1520
1521 // Create Rs functions.
1522 rsFunc_.reserve(rec.size());
1523
1524 auto getArray = [](const std::vector<double>& input)
1525 {
1526 if constexpr (std::is_same_v<Scalar,double>) {
1527 return input;
1528 } else {
1529 std::vector<Scalar> output;
1530 output.resize(input.size());
1531 std::ranges::copy(input, output.begin());
1532 return output;
1533 }
1534 };
1535
1536 if (FluidSystem::enableDissolvedGas()) {
1537 for (std::size_t i = 0; i < rec.size(); ++i) {
1538 if (eqlmap.cells(i).empty()) {
1539 rsFunc_.push_back(std::shared_ptr<Miscibility::RsVD<FluidSystem>>());
1540 continue;
1541 }
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());
1550 rsFunc_.push_back(std::make_shared<Miscibility::RsVD<FluidSystem>>(pvtIdx,
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());
1556 rsFunc_.push_back(std::make_shared<Miscibility::PBVD<FluidSystem>>(pvtIdx,
1557 depthColumn, pbubColumn));
1558
1559 } else {
1560 throw std::runtime_error("Cannot initialise: RSVD or PBVD table not available.");
1561 }
1562
1563 }
1564 else {
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.");
1569 }
1570 const Scalar pContact = rec[i].datumDepthPressure();
1571 const Scalar TContact = 273.15 + 20; // standard temperature for now
1572 rsFunc_.push_back(std::make_shared<Miscibility::RsSatAtContact<FluidSystem>>(pvtIdx, pContact, TContact));
1573 }
1574 }
1575 }
1576 else {
1577 for (std::size_t i = 0; i < rec.size(); ++i) {
1578 rsFunc_.push_back(std::make_shared<Miscibility::NoMixing<Scalar>>());
1579 }
1580 }
1581
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()) {
1586 rvFunc_.push_back(std::shared_ptr<Miscibility::RvVD<FluidSystem>>());
1587 continue;
1588 }
1589 const int pvtIdx = regionPvtIdx_[i];
1590 if (!rec[i].wetGasInitConstantRv()) {
1591 const TableContainer& rvvdTables = tables.getRvvdTables();
1592 const TableContainer& pdvdTables = tables.getPdvdTables();
1593
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());
1598 rvFunc_.push_back(std::make_shared<Miscibility::RvVD<FluidSystem>>(pvtIdx,
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());
1604 rvFunc_.push_back(std::make_shared<Miscibility::PDVD<FluidSystem>>(pvtIdx,
1605 depthColumn, pdewColumn));
1606 } else {
1607 throw std::runtime_error("Cannot initialise: RVVD or PDCD table not available.");
1608 }
1609 }
1610 else {
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.");
1616 }
1617 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1618 const Scalar TContact = 273.15 + 20; // standard temperature for now
1619 rvFunc_.push_back(std::make_shared<Miscibility::RvSatAtContact<FluidSystem>>(pvtIdx,pContact, TContact));
1620 }
1621 }
1622 }
1623 else {
1624 for (std::size_t i = 0; i < rec.size(); ++i) {
1625 rvFunc_.push_back(std::make_shared<Miscibility::NoMixing<Scalar>>());
1626 }
1627 }
1628
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()) {
1633 rvwFunc_.push_back(std::shared_ptr<Miscibility::RvwVD<FluidSystem>>());
1634 continue;
1635 }
1636 const int pvtIdx = regionPvtIdx_[i];
1637 if (!rec[i].humidGasInitConstantRvw()) {
1638 const TableContainer& rvwvdTables = tables.getRvwvdTables();
1639
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());
1644 rvwFunc_.push_back(std::make_shared<Miscibility::RvwVD<FluidSystem>>(pvtIdx,
1645 depthColumn, rvwvdColumn));
1646 } else {
1647 throw std::runtime_error("Cannot initialise: RVWVD table not available.");
1648 }
1649 }
1650 else {
1651 const auto oilActive = FluidSystem::phaseIsActive(FluidSystem::oilPhaseIdx);
1652 if (oilActive) {
1653 if (rec[i].gasOilContactDepth() != rec[i].datumDepth()) {
1654 rvwFunc_.push_back(std::make_shared<Miscibility::NoMixing<Scalar>>());
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);
1659 } else {
1660 // pg = po + Pcgo = po + (pg - po)
1661 // for gas-condensate with initial no oil zone: water-oil contact depth (OWC) equal gas-oil contact depth (GOC)
1662 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].gasOilContactCapillaryPressure();
1663 const Scalar TContact = 273.15 + 20; // standard temperature for now
1664 rvwFunc_.push_back(std::make_shared<Miscibility::RvwSatAtContact<FluidSystem>>(pvtIdx,pContact, TContact));
1665 }
1666 }
1667 else {
1668 // two-phase gas-water sytem: water-oil contact depth is taken equal to gas-water contact depth (GWC)
1669 // and water-oil capillary pressure (Pcwo) is taken equal to gas-water capillary pressure (Pcgw) at GWC
1670 if (rec[i].waterOilContactDepth() != rec[i].datumDepth()) {
1671 rvwFunc_.push_back(std::make_shared<Miscibility::NoMixing<Scalar>>());
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);
1676 } else {
1677 // pg = pw + Pcgw = pw + (pg - pw)
1678 const Scalar pContact = rec[i].datumDepthPressure() + rec[i].waterOilContactCapillaryPressure();
1679 const Scalar TContact = 273.15 + 20; // standard temperature for now
1680 rvwFunc_.push_back(std::make_shared<Miscibility::RvwSatAtContact<FluidSystem>>(pvtIdx,pContact, TContact));
1681 }
1682 }
1683 }
1684 }
1685 }
1686 else {
1687 for (std::size_t i = 0; i < rec.size(); ++i) {
1688 rvwFunc_.push_back(std::make_shared<Miscibility::NoMixing<Scalar>>());
1689 }
1690 }
1691
1692
1693 // EXTRACT the initial temperature
1694 updateInitialTemperature_(eclipseState, eqlmap);
1695
1696 // EXTRACT the initial salt concentration
1697 updateInitialSaltConcentration_(eclipseState, eqlmap);
1698
1699 // EXTRACT the initial salt saturation
1700 updateInitialSaltSaturation_(eclipseState, eqlmap);
1701
1702 // Compute pressures, saturations, rs and rv factors.
1703 const auto& comm = grid.comm();
1704 calcPressSatRsRv(eqlmap, rec, materialLawManager, gridView, comm, grav);
1705
1706 // modify the pressure and saturation for numerical aquifer cells
1707 applyNumericalAquifers_(gridView, num_aquifers, eclipseState.runspec().co2Storage() || eclipseState.runspec().h2Storage());
1708
1709 // Modify oil pressure in no-oil regions so that the pressures of present phases can
1710 // be recovered from the oil pressure and capillary relations.
1711}
1712
1713template<class FluidSystem,
1714 class Grid,
1715 class GridView,
1716 class ElementMapper,
1717 class CartesianIndexMapper>
1718template<class RMap>
1719void InitialStateComputer<FluidSystem,
1720 Grid,
1721 GridView,
1722 ElementMapper,
1723 CartesianIndexMapper>::
1724updateInitialTemperature_(const EclipseState& eclState, const RMap& reg)
1725{
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);
1735 }
1736 } else {
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, /*extrapolate=*/true);
1745 }
1746 }
1747 }
1748}
1749
1750template<class FluidSystem,
1751 class Grid,
1752 class GridView,
1753 class ElementMapper,
1754 class CartesianIndexMapper>
1755template<class RMap>
1756void InitialStateComputer<FluidSystem,
1757 Grid,
1758 GridView,
1759 ElementMapper,
1760 CartesianIndexMapper>::
1761updateInitialSaltConcentration_(const EclipseState& eclState, const RMap& reg)
1762{
1763 const int numEquilReg = rsFunc_.size();
1764 saltVdTable_.resize(numEquilReg);
1765 const auto& tables = eclState.getTableManager();
1766 const TableContainer& saltvdTables = tables.getSaltvdTables();
1767
1768 // If no saltvd table is given, we create a trivial table for the density calculations
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);
1774 }
1775 } else {
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());
1779
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, /*extrapolate=*/true);
1784 }
1785 }
1786 }
1787}
1788
1789template<class FluidSystem,
1790 class Grid,
1791 class GridView,
1792 class ElementMapper,
1793 class CartesianIndexMapper>
1794template<class RMap>
1795void InitialStateComputer<FluidSystem,
1796 Grid,
1797 GridView,
1798 ElementMapper,
1799 CartesianIndexMapper>::
1800updateInitialSaltSaturation_(const EclipseState& eclState, const RMap& reg)
1801{
1802 const int numEquilReg = rsFunc_.size();
1803 saltpVdTable_.resize(numEquilReg);
1804 const auto& tables = eclState.getTableManager();
1805 const TableContainer& saltpvdTables = tables.getSaltpvdTables();
1806
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());
1810
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, /*extrapolate=*/true);
1815 }
1816 }
1817}
1818
1819template<class FluidSystem,
1820 class Grid,
1821 class GridView,
1822 class ElementMapper,
1823 class CartesianIndexMapper>
1824void InitialStateComputer<FluidSystem,
1825 Grid,
1826 GridView,
1827 ElementMapper,
1828 CartesianIndexMapper>::
1829updateCellProps_(const GridView& gridView,
1830 const NumericalAquifers& aquifer)
1831{
1832 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1833 int numElements = gridView.size(/*codim=*/0);
1834 cellCenterDepth_.resize(numElements);
1835 cellCenterXY_.resize(numElements);
1836 cellCorners_.resize(numElements);
1837 cellZSpan_.resize(numElements);
1838 cellZMinMax_.resize(numElements);
1839
1840 auto elemIt = gridView.template begin</*codim=*/0>();
1841 const auto& elemEndIt = gridView.template end</*codim=*/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;
1862 }
1863 }
1864 }
1865}
1866
1867template<class FluidSystem,
1868 class Grid,
1869 class GridView,
1870 class ElementMapper,
1871 class CartesianIndexMapper>
1872void InitialStateComputer<FluidSystem,
1873 Grid,
1874 GridView,
1875 ElementMapper,
1876 CartesianIndexMapper>::
1877applyNumericalAquifers_(const GridView& gridView,
1878 const NumericalAquifers& aquifer,
1879 const bool co2store_or_h2store)
1880{
1881 const auto num_aqu_cells = aquifer.allAquiferCells();
1882 if (num_aqu_cells.empty()) return;
1883
1884 // Check if water phase is active, or in the case of CO2STORE and H2STORE, water is modelled as oil phase
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" };
1889 }
1890
1891 ElementMapper elemMapper(gridView, Dune::mcmgElementLayout());
1892 auto elemIt = gridView.template begin</*codim=*/0>();
1893 const auto& elemEndIt = gridView.template end</*codim=*/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()) {
1902 // numerical aquifer cells are filled with water initially
1903 this->sat_[watPos][elemIdx] = 1.;
1904
1905 if (!co2store_or_h2store && FluidSystem::phaseIsActive(oilPos)) {
1906 this->sat_[oilPos][elemIdx] = 0.;
1907 }
1908
1909 if (FluidSystem::phaseIsActive(gasPos)) {
1910 this->sat_[gasPos][elemIdx] = 0.;
1911 }
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);
1916 OpmLog::info(msg);
1917
1918 // if pressure is specified for numerical aquifers, we use these pressure values
1919 // for numerical aquifer cells
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;
1925 }
1926 if (FluidSystem::phaseIsActive(oilPos)) {
1927 this->pp_[oilPos][elemIdx] = pres;
1928 }
1929 }
1930 }
1931 }
1932}
1933
1934template<class FluidSystem,
1935 class Grid,
1936 class GridView,
1937 class ElementMapper,
1938 class CartesianIndexMapper>
1939template<class RMap>
1940void InitialStateComputer<FluidSystem,
1941 Grid,
1942 GridView,
1943 ElementMapper,
1944 CartesianIndexMapper>::
1945setRegionPvtIdx(const EclipseState& eclState, const RMap& reg)
1946{
1947 const auto& pvtnumData = eclState.fieldProps().get_int("PVTNUM");
1948
1949 for (const auto& r : reg.activeRegions()) {
1950 const auto& cells = reg.cells(r);
1951 regionPvtIdx_[r] = pvtnumData[*cells.begin()] - 1;
1952 }
1953}
1954
1955template<class FluidSystem,
1956 class Grid,
1957 class GridView,
1958 class ElementMapper,
1959 class CartesianIndexMapper>
1960template<class RMap, class MaterialLawManager, class Comm>
1961void InitialStateComputer<FluidSystem,
1962 Grid,
1963 GridView,
1964 ElementMapper,
1965 CartesianIndexMapper>::
1966calcPressSatRsRv(const RMap& reg,
1967 const std::vector<EquilRecord>& rec,
1968 MaterialLawManager& materialLawManager,
1969 const GridView& gridView,
1970 const Comm& comm,
1971 const Scalar grav)
1972{
1973 using PhaseSat = Details::PhaseSaturations<
1974 MaterialLawManager, FluidSystem, EquilReg<Scalar>, typename RMap::CellId
1975 >;
1976
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>{};
1980
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);
1984
1985 Details::verticalExtent(cells, cellZMinMax_, comm, vspan);
1986
1987 const auto acc = rec[r].initializationTargetAccuracy();
1988 if (acc > 0) {
1989 // The grid blocks are treated as being tilted
1990 // First check if the region has cells
1991 if (cells.empty()) {
1992 regionIsEmpty[r] = 1;
1993 continue;
1994 }
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]
1998 };
1999 // Ensure contacts are within the span
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);
2003 // For titled blocks, we can use a simple weightening based on title of the grid
2004 // this->equilibrateTiltedFaultBlockSimple(cells, eqreg, gridView, acc, ptable, psat);
2005 this->equilibrateTiltedFaultBlock(cells, eqreg, gridView, acc, ptable, psat);
2006 }
2007 else if (acc == 0) {
2008 if (cells.empty()) {
2009 regionIsEmpty[r] = 1;
2010 continue;
2011 }
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]
2015 };
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);
2019 // Centre-point method
2020 this->equilibrateCellCentres(cells, eqreg, ptable, psat);
2021 }
2022 else if (acc < 0) {
2023 if (cells.empty()) {
2024 regionIsEmpty[r] = 1;
2025 continue;
2026 }
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]
2030 };
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);
2034 // Horizontal subdivision
2035 this->equilibrateHorizontal(cells, eqreg, -acc, ptable, psat);
2036 }
2037 }
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]) //region is empty on all partitions
2042 OpmLog::warning("Equilibration region " + std::to_string(r + 1)
2043 + " has no active cells");
2044 }
2045 }
2046}
2047
2048template<class FluidSystem,
2049 class Grid,
2050 class GridView,
2051 class ElementMapper,
2052 class CartesianIndexMapper>
2053template<class CellRange, class EquilibrationMethod>
2054void InitialStateComputer<FluidSystem,
2055 Grid,
2056 GridView,
2057 ElementMapper,
2058 CartesianIndexMapper>::
2059cellLoop(const CellRange& cells,
2060 EquilibrationMethod&& eqmethod)
2061{
2062 const auto oilPos = FluidSystem::oilPhaseIdx;
2063 const auto gasPos = FluidSystem::gasPhaseIdx;
2064 const auto watPos = FluidSystem::waterPhaseIdx;
2065
2066 const auto oilActive = FluidSystem::phaseIsActive(oilPos);
2067 const auto gasActive = FluidSystem::phaseIsActive(gasPos);
2068 const auto watActive = FluidSystem::phaseIsActive(watPos);
2069
2070 auto pressures = Details::PhaseQuantityValue<Scalar>{};
2071 auto saturations = Details::PhaseQuantityValue<Scalar>{};
2072 Scalar Rs = 0.0;
2073 Scalar Rv = 0.0;
2074 Scalar Rvw = 0.0;
2075
2076 for (const auto& cell : cells) {
2077 eqmethod(cell, pressures, saturations, Rs, Rv, Rvw);
2078
2079 if (oilActive) {
2080 this->pp_ [oilPos][cell] = pressures.oil;
2081 this->sat_[oilPos][cell] = saturations.oil;
2082 }
2083
2084 if (gasActive) {
2085 this->pp_ [gasPos][cell] = pressures.gas;
2086 this->sat_[gasPos][cell] = saturations.gas;
2087 }
2088
2089 if (watActive) {
2090 this->pp_ [watPos][cell] = pressures.water;
2091 this->sat_[watPos][cell] = saturations.water;
2092 }
2093
2094 if (oilActive && gasActive) {
2095 this->rs_[cell] = Rs;
2096 this->rv_[cell] = Rv;
2097 }
2098
2099 if (watActive && gasActive) {
2100 this->rvw_[cell] = Rvw;
2101 }
2102 }
2103}
2104
2105template<class FluidSystem,
2106 class Grid,
2107 class GridView,
2108 class ElementMapper,
2109 class CartesianIndexMapper>
2110template<class CellRange, class PressTable, class PhaseSat>
2111void InitialStateComputer<FluidSystem,
2112 Grid,
2113 GridView,
2114 ElementMapper,
2115 CartesianIndexMapper>::
2116equilibrateCellCentres(const CellRange& cells,
2117 const EquilReg<Scalar>& eqreg,
2118 const PressTable& ptable,
2119 PhaseSat& psat)
2120{
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]
2125 (const CellID cell,
2126 Details::PhaseQuantityValue<Scalar>& pressures,
2127 Details::PhaseQuantityValue<Scalar>& saturations,
2128 Scalar& Rs,
2129 Scalar& Rv,
2130 Scalar& Rvw) -> void
2131 {
2132 const auto pos = CellPos {
2133 cell, cellCenterDepth_[cell]
2134 };
2135
2136 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2137 pressures = psat.correctedPhasePressures();
2138
2139 const auto temp = this->temperature_[cell];
2140
2141 Rs = eqreg.dissolutionCalculator()
2142 (pos.depth, pressures.oil, temp, saturations.gas);
2143
2144 Rv = eqreg.evaporationCalculator()
2145 (pos.depth, pressures.gas, temp, saturations.oil);
2146
2147 Rvw = eqreg.waterEvaporationCalculator()
2148 (pos.depth, pressures.gas, temp, saturations.water);
2149 });
2150}
2151
2152template<class FluidSystem,
2153 class Grid,
2154 class GridView,
2155 class ElementMapper,
2156 class CartesianIndexMapper>
2157template<class CellRange, class PressTable, class PhaseSat>
2158void InitialStateComputer<FluidSystem,
2159 Grid,
2160 GridView,
2161 ElementMapper,
2162 CartesianIndexMapper>::
2163equilibrateHorizontal(const CellRange& cells,
2164 const EquilReg<Scalar>& eqreg,
2165 const int acc,
2166 const PressTable& ptable,
2167 PhaseSat& psat)
2168{
2169 using CellPos = typename PhaseSat::Position;
2170 using CellID = std::remove_cv_t<std::remove_reference_t<
2171 decltype(std::declval<CellPos>().cell)>>;
2172
2173 this->cellLoop(cells, [this, acc, &eqreg, &ptable, &psat]
2174 (const CellID cell,
2175 Details::PhaseQuantityValue<Scalar>& pressures,
2176 Details::PhaseQuantityValue<Scalar>& saturations,
2177 Scalar& Rs,
2178 Scalar& Rv,
2179 Scalar& Rvw) -> void
2180 {
2181 pressures .reset();
2182 saturations.reset();
2183
2184 Scalar totfrac = 0.0;
2185 for (const auto& [depth, frac] : Details::horizontalSubdivision(cell, cellZSpan_[cell], acc)) {
2186 const auto pos = CellPos { cell, depth };
2187
2188 saturations.axpy(psat.deriveSaturations(pos, eqreg, ptable), frac);
2189 pressures .axpy(psat.correctedPhasePressures(), frac);
2190
2191 totfrac += frac;
2192 }
2193
2194 if (totfrac > 0.) {
2195 saturations /= totfrac;
2196 pressures /= totfrac;
2197 } else {
2198 // Fall back to centre point method for zero-thickness cells.
2199 const auto pos = CellPos {
2200 cell, cellCenterDepth_[cell]
2201 };
2202
2203 saturations = psat.deriveSaturations(pos, eqreg, ptable);
2204 pressures = psat.correctedPhasePressures();
2205 }
2206
2207 const auto temp = this->temperature_[cell];
2208 const auto cz = cellCenterDepth_[cell];
2209
2210 Rs = eqreg.dissolutionCalculator()
2211 (cz, pressures.oil, temp, saturations.gas);
2212
2213 Rv = eqreg.evaporationCalculator()
2214 (cz, pressures.gas, temp, saturations.oil);
2215
2216 Rvw = eqreg.waterEvaporationCalculator()
2217 (cz, pressures.gas, temp, saturations.water);
2218 });
2219}
2220
2221template<class FluidSystem, class Grid, class GridView, 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,
2227 const int acc,
2228 const PressTable& ptable,
2229 PhaseSat& psat)
2230{
2231 using CellPos = typename PhaseSat::Position;
2232 using CellID = std::remove_cv_t<std::remove_reference_t<
2233 decltype(std::declval<CellPos>().cell)>>;
2234
2235 this->cellLoop(cells, [this, acc, &eqreg, &ptable, &psat, &gridView]
2236 (const CellID cell,
2237 Details::PhaseQuantityValue<Scalar>& pressures,
2238 Details::PhaseQuantityValue<Scalar>& saturations,
2239 Scalar& Rs,
2240 Scalar& Rv,
2241 Scalar& Rvw) -> void
2242 {
2243 pressures.reset();
2244 saturations.reset();
2245 Scalar totalWeight = 0.0;
2246
2247 // We assume grid blocks are treated as being tilted
2248 const auto& [zmin, zmax] = cellZMinMax_[cell];
2249 const Scalar cellThickness = zmax - zmin;
2250 const Scalar halfThickness = cellThickness / 2.0;
2251
2252 // Calculate dip parameters from corner point geometry
2253 Scalar dipAngle, dipAzimuth;
2254 Details::computeBlockDip(cellCorners_[cell], dipAngle, dipAzimuth);
2255
2256 // Reference point for TVD calculations
2257 std::array<Scalar, 3> referencePoint = {
2258 cellCenterXY_[cell].first,
2259 cellCenterXY_[cell].second,
2260 cellCenterDepth_[cell]
2261 };
2262
2263 // We have acc levels within each half (upper and lower) of the block
2264 const int numLevelsPerHalf = std::min(20, acc);
2265
2266 // Create subdivisions for upper and lower halves with cross-section weighting
2267 std::vector<std::pair<Scalar, Scalar>> levels;
2268
2269 // Subdivide upper and lower halves separately
2270 for (int side = 0; side < 2; ++side) {
2271 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2272
2273 for (int i = 0; i < numLevelsPerHalf; ++i) {
2274 // Calculate depth at the center of this subdivision
2275 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2276
2277 // A simple way: we can estimate cross-section weight based on dip angle
2278 // For horizontal cells: weight = 1.0, for tilted cells: weight decreases with dip
2279 Scalar crossSectionWeight = (halfThickness / numLevelsPerHalf);
2280
2281 // Apply dip correction to weight (cross-section area decreases with dip)
2282 if (std::abs(dipAngle) > 1e-10) {
2283 crossSectionWeight /= std::cos(dipAngle);
2284 }
2285
2286 levels.emplace_back(depth, crossSectionWeight);
2287 }
2288 }
2289
2290 for (const auto& [depth, weight] : levels) {
2291 // Convert measured depth to True Vertical Depth for tilted blocks
2292 const auto& [x, y] = cellCenterXY_[cell];
2294 depth, x, y, dipAngle, dipAzimuth, referencePoint);
2295
2296 const auto pos = CellPos{cell, tvd};
2297
2298 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2299 auto localPressures = psat.correctedPhasePressures();
2300
2301 // Apply cross-section weighted averaging
2302 saturations.axpy(localSaturations, weight);
2303 pressures.axpy(localPressures, weight);
2304 totalWeight += weight;
2305 }
2306
2307 // Normalize results
2308 if (totalWeight > 1e-10) {
2309 saturations /= totalWeight;
2310 pressures /= totalWeight;
2311 } else {
2312 // Fallback to center point method using TVD
2313 const auto& [x, y] = cellCenterXY_[cell];
2314 Scalar tvdCenter = Details::calculateTrueVerticalDepth(
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();
2319 }
2320
2321 // Compute solution ratios at cell center TVD
2322 const auto temp = this->temperature_[cell];
2323 const auto& [x, y] = cellCenterXY_[cell];
2324 Scalar tvdCenter = Details::calculateTrueVerticalDepth(
2325 cellCenterDepth_[cell], x, y, dipAngle, dipAzimuth, referencePoint);
2326
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);
2330 });
2331}
2332
2333template<class FluidSystem, class Grid, class GridView, 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,
2339 const int acc,
2340 const PressTable& ptable,
2341 PhaseSat& psat)
2342{
2343 using CellPos = typename PhaseSat::Position;
2344 using CellID = std::remove_cv_t<std::remove_reference_t<
2345 decltype(std::declval<CellPos>().cell)>>;
2346
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;
2351 }
2352
2353 // Face Area Calculation
2354 auto polygonArea = [](const std::vector<std::array<Scalar, 2>>& pts) {
2355 if (pts.size() < 3) return Scalar(0);
2356 Scalar area = 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];
2360 }
2361 return std::abs(area) * Scalar(0.5);
2362 };
2363
2364 // Compute horizontal cross-section at given depth
2365 auto computeCrossSectionArea = [&](const CellID cell, Scalar depth) -> Scalar {
2366 try {
2367 const auto& entity = entityMap[cell];
2368 const auto& geometry = entity.geometry();
2369 const int numCorners = geometry.corners();
2370
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])};
2375 }
2376
2377 // Find all intersections between horizontal plane and cell edges
2378 std::vector<std::array<Scalar, 2>> intersectionPoints;
2379 const Scalar tol = 1e-10;
2380
2381 // Check all edges between corners (could be optimized further)
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];
2386
2387 if ((za - depth) * (zb - depth) <= 0.0 && std::abs(za - zb) > tol) {
2388 // Edge crosses the horizontal plane
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});
2393 }
2394 }
2395 }
2396
2397 // Remove duplicates
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;
2401 };
2402
2403 intersectionPoints.erase(
2404 std::unique(intersectionPoints.begin(), intersectionPoints.end(), pointsEqual),
2405 intersectionPoints.end()
2406 );
2407 }
2408
2409 if (intersectionPoints.size() < 3) {
2410 // No valid grid found, use fallback
2411 return 0.0;
2412 }
2413
2414 // Order points counter-clockwise around centroid
2415 Scalar cx = 0, cy = 0;
2416 for (const auto& p : intersectionPoints) {
2417 cx += p[0]; cy += p[1];
2418 }
2419 cx /= intersectionPoints.size();
2420 cy /= intersectionPoints.size();
2421
2422 // Sorting
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);
2425 };
2426
2427 std::ranges::sort(intersectionPoints, angleCompare);
2428
2429 return polygonArea(intersectionPoints);
2430
2431 } catch (const std::exception& e) {
2432 return 0.0;
2433 }
2434 };
2435
2436 auto cellProcessor = [this, acc, &eqreg, &ptable, &psat, &computeCrossSectionArea]
2437 (const CellID cell,
2438 Details::PhaseQuantityValue<Scalar>& pressures,
2439 Details::PhaseQuantityValue<Scalar>& saturations,
2440 Scalar& Rs,
2441 Scalar& Rv,
2442 Scalar& Rvw) -> void
2443 {
2444 pressures.reset();
2445 saturations.reset();
2446 Scalar totalWeight = 0.0;
2447
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;
2452
2453 // Calculate dip parameters from corner point geometry
2454 Scalar dipAngle, dipAzimuth;
2455 Details::computeBlockDip(this->cellCorners_[cell], dipAngle, dipAzimuth);
2456
2457 // Reference point for TVD calculations
2458 std::array<Scalar, 3> referencePoint = {
2459 this->cellCenterXY_[cell].first,
2460 this->cellCenterXY_[cell].second,
2461 cellCenterDepth_[cell]
2462 };
2463
2464 // We have acc levels within each half (upper and lower) of the block
2465 const int numLevelsPerHalf = std::min(20, acc);
2466
2467 // Create subdivisions for upper and lower halves with cross-section weighting
2468 std::vector<std::pair<Scalar, Scalar>> levels;
2469
2470 // Subdivide upper and lower halves separately
2471 for (int side = 0; side < 2; ++side) {
2472 Scalar halfStart = (side == 0) ? zmin : zmin + halfThickness;
2473
2474 for (int i = 0; i < numLevelsPerHalf; ++i) {
2475 // Calculate depth at the center of this subdivision
2476 Scalar depth = halfStart + (i + 0.5) * (halfThickness / numLevelsPerHalf);
2477
2478 // Compute cross-section area at this depth
2479 Scalar crossSectionArea = computeCrossSectionArea(cell, depth);
2480
2481 // Weight is proportional to: area × Δz (volume element)
2482 Scalar weight = crossSectionArea * (halfThickness / numLevelsPerHalf);
2483
2484 levels.emplace_back(depth, weight);
2485 }
2486 }
2487
2488 // Maybe not necessary (for debug)
2489 bool hasValidAreas = false;
2490 for (const auto& level : levels) {
2491 if (level.second > 1e-10) {
2492 hasValidAreas = true;
2493 break;
2494 }
2495 }
2496
2497 if (!hasValidAreas) {
2498 // Fallback to dip-based weighting as used in equilibrateTiltedFaultBlockSimple
2499 levels.clear();
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);
2507 }
2508 levels.emplace_back(depth, weight);
2509 }
2510 }
2511 }
2512
2513 for (const auto& level : levels) {
2514 Scalar depth = level.first;
2515 Scalar weight = level.second;
2516
2517 // Convert measured depth to True Vertical Depth for tilted blocks
2518 const auto& xy = this->cellCenterXY_[cell];
2520 depth, xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2521
2522 const auto pos = CellPos{cell, tvd};
2523
2524 auto localSaturations = psat.deriveSaturations(pos, eqreg, ptable);
2525 auto localPressures = psat.correctedPhasePressures();
2526
2527 // Apply cross-section weighted averaging
2528 saturations.axpy(localSaturations, weight);
2529 pressures.axpy(localPressures, weight);
2530 totalWeight += weight;
2531 }
2532
2533 if (totalWeight > 1e-10) {
2534 saturations /= totalWeight;
2535 pressures /= totalWeight;
2536 } else {
2537 // Fallback to center point method using TVD
2538 const auto& xy = this->cellCenterXY_[cell];
2539 Scalar tvdCenter = Details::calculateTrueVerticalDepth(
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();
2544 }
2545
2546 // Compute solution ratios at cell center TVD
2547 const auto temp = this->temperature_[cell];
2548 const auto& xy = this->cellCenterXY_[cell];
2549 Scalar tvdCenter = Details::calculateTrueVerticalDepth(
2550 this->cellCenterDepth_[cell], xy.first, xy.second, dipAngle, dipAzimuth, referencePoint);
2551
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);
2555 };
2556
2557 this->cellLoop(cells, cellProcessor);
2558}
2559}
2560} // namespace EQUIL
2561} // namespace Opm
2562
2563#endif // OPM_INIT_STATE_EQUIL_IMPL_HPP
#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 &reg, 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 &reg, 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:439
Definition: EquilibrationHelpers.hpp:162
Definition: EquilibrationHelpers.hpp:500
Definition: EquilibrationHelpers.hpp:322
Definition: EquilibrationHelpers.hpp:560
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