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