FluidMatrixInteractionBlackoil.hpp
Go to the documentation of this file.
1/*
2 Copyright 2010 SINTEF ICT, Applied Mathematics.
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
20#ifndef OPM_FLUIDMATRIXINTERACTIONBLACKOIL_HEADER_INCLUDED
21#define OPM_FLUIDMATRIXINTERACTIONBLACKOIL_HEADER_INCLUDED
22
23#include <opm/core/utility/UniformTableLinear.hpp>
24#include <opm/core/utility/buildUniformMonotoneTable.hpp>
25#include <opm/parser/eclipse/Deck/Deck.hpp>
26#include <opm/parser/eclipse/Parser/ParseContext.hpp>
27#include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
28#include <opm/parser/eclipse/EclipseState/Tables/SwofTable.hpp>
29#include <opm/parser/eclipse/EclipseState/Tables/SgofTable.hpp>
30#include <opm/parser/eclipse/EclipseState/Tables/TableManager.hpp>
31#include "BlackoilDefs.hpp"
32#include <iostream>
33#include <stdexcept>
34
35namespace Opm
36{
37
38// Forward declaration needed by associated parameters class.
39template <class ScalarT, class ParamsT>
40class FluidMatrixInteractionBlackoil;
41
42template <class ScalarT>
44{
45public:
46 typedef ScalarT Scalar;
47 void init(Opm::DeckConstPtr deck)
48 {
49 ParseContext parseContext;
50 EclipseState eclipseState(deck , parseContext);
51 // Extract input data.
52 const auto& tables = eclipseState.getTableManager();
53 const auto& swofTable = tables->getSwofTables().getTable<SwofTable>(0);
54 const auto& sgofTable = tables->getSgofTables().getTable<SgofTable>(0);
55
56 std::vector<double> sw = swofTable.getSwColumn().vectorCopy();
57 std::vector<double> krw = swofTable.getKrwColumn().vectorCopy();
58 std::vector<double> krow = swofTable.getKrowColumn().vectorCopy();
59 std::vector<double> pcow = swofTable.getPcowColumn().vectorCopy();
60 std::vector<double> sg = sgofTable.getSgColumn().vectorCopy();
61 std::vector<double> krg = sgofTable.getKrgColumn().vectorCopy();
62 std::vector<double> krog = sgofTable.getKrogColumn().vectorCopy();
63 std::vector<double> pcog = sgofTable.getPcogColumn().vectorCopy();
64
65 // Create tables for krw, krow, krg and krog.
66 int samples = 200;
67 buildUniformMonotoneTable(sw, krw, samples, krw_);
68 buildUniformMonotoneTable(sw, krow, samples, krow_);
69 buildUniformMonotoneTable(sg, krg, samples, krg_);
70 buildUniformMonotoneTable(sg, krog, samples, krog_);
71 krocw_ = krow[0]; // At connate water -> ecl. SWOF
72
73 // Create tables for pcow and pcog.
74 buildUniformMonotoneTable(sw, pcow, samples, pcow_);
75 buildUniformMonotoneTable(sg, pcog, samples, pcog_);
76 }
77
78private:
79 template <class S, class P>
81
82 Opm::utils::UniformTableLinear<Scalar> krw_;
83 Opm::utils::UniformTableLinear<Scalar> krow_;
84 Opm::utils::UniformTableLinear<Scalar> pcow_;
85 Opm::utils::UniformTableLinear<Scalar> krg_;
86 Opm::utils::UniformTableLinear<Scalar> krog_;
87 Opm::utils::UniformTableLinear<Scalar> pcog_;
88 Scalar krocw_; // = krow_(s_wc)
89};
90
91
97template <class ScalarT, class ParamsT = FluidMatrixInteractionBlackoilParams<ScalarT> >
99{
100public:
101 typedef ParamsT Params;
102 typedef typename Params::Scalar Scalar;
103
114 template <class pcContainerT, class SatContainerT>
115 static void pC(pcContainerT &pc,
116 const Params &params,
117 const SatContainerT &saturations,
118 Scalar /*temperature*/)
119 {
120 Scalar sw = saturations[Aqua];
121 Scalar sg = saturations[Vapour];
122 pc[Liquid] = 0.0;
123 pc[Aqua] = params.pcow_(sw);
124 pc[Vapour] = params.pcog_(sg);
125 }
126
138 template <class SatContainerT, class pcContainerT>
139 static void S(SatContainerT &saturations,
140 const Params &params,
141 const pcContainerT &pc,
142 Scalar /*temperature*/)
143 {
144 std::cerr << "FluidMatrixInteractionBlackoil::S() is not implemented yet\n";
145 throw std::logic_error("Not implemented");
146 }
147
148
152 template <class krContainerT, class SatContainerT>
153 static void kr(krContainerT &kr,
154 const Params &params,
155 const SatContainerT &saturations,
156 Scalar /*temperature*/)
157 {
158 // Stone-II relative permeability model.
159 Scalar sw = saturations[Aqua];
160 Scalar sg = saturations[Vapour];
161 Scalar krw = params.krw_(sw);
162 Scalar krg = params.krg_(sg);
163 Scalar krow = params.krow_(sw);
164 Scalar krog = params.krog_(sg);
165 Scalar krocw = params.krocw_;
166 kr[Aqua] = krw;
167 kr[Vapour] = krg;
168 kr[Liquid] = krocw*((krow/krocw + krw)*(krog/krocw + krg) - krw - krg);
169 if (kr[Liquid] < 0.0) {
170 kr[Liquid] = 0.0;
171 }
172 }
173
174
178 template <class krContainerT, class SatContainerT>
179 static void dkr(krContainerT &dkr,
180 const Params &params,
181 const SatContainerT &saturations,
182 Scalar /*temperature*/)
183 {
184 for (int p1 = 0; p1 < numPhases; ++p1) {
185 for (int p2 = 0; p2 < numPhases; ++p2) {
186 dkr[p1][p2] = 0.0;
187 }
188 }
189 // Stone-II relative permeability model.
190 Scalar sw = saturations[Aqua];
191 Scalar sg = saturations[Vapour];
192 Scalar krw = params.krw_(sw);
193 Scalar dkrww = params.krw_.derivative(sw);
194 Scalar krg = params.krg_(sg);
195 Scalar dkrgg = params.krg_.derivative(sg);
196 Scalar krow = params.krow_(sw);
197 Scalar dkrow = params.krow_.derivative(sw);
198 Scalar krog = params.krog_(sg);
199 Scalar dkrog = params.krog_.derivative(sg);
200 Scalar krocw = params.krocw_;
201 dkr[Aqua][Aqua] = dkrww;
202 dkr[Vapour][Vapour] = dkrgg;
203 dkr[Liquid][Aqua] = krocw*((dkrow/krocw + dkrww)*(krog/krocw + krg) - dkrww);
204 dkr[Liquid][Vapour] = krocw*((krow/krocw + krw)*(dkrog/krocw + dkrgg) - dkrgg);
205 }
206};
207
208} // namespace Opm
209
210
211
212
213#endif // OPM_FLUIDMATRIXINTERACTIONBLACKOIL_HEADER_INCLUDED
Definition: BlackoilDefs.hpp:31
double Scalar
Definition: BlackoilDefs.hpp:39
@ Aqua
Definition: BlackoilDefs.hpp:37
@ Vapour
Definition: BlackoilDefs.hpp:37
@ Liquid
Definition: BlackoilDefs.hpp:37
@ numPhases
Definition: BlackoilDefs.hpp:34
Definition: FluidMatrixInteractionBlackoil.hpp:44
void init(Opm::DeckConstPtr deck)
Definition: FluidMatrixInteractionBlackoil.hpp:47
ScalarT Scalar
Definition: FluidMatrixInteractionBlackoil.hpp:46
Capillary pressures and relative permeabilities for a black oil system.
Definition: FluidMatrixInteractionBlackoil.hpp:99
static void kr(krContainerT &kr, const Params &params, const SatContainerT &saturations, Scalar)
The relative permeability of all phases.
Definition: FluidMatrixInteractionBlackoil.hpp:153
static void dkr(krContainerT &dkr, const Params &params, const SatContainerT &saturations, Scalar)
The saturation derivatives of relative permeability of all phases.
Definition: FluidMatrixInteractionBlackoil.hpp:179
Params::Scalar Scalar
Definition: FluidMatrixInteractionBlackoil.hpp:102
static void pC(pcContainerT &pc, const Params &params, const SatContainerT &saturations, Scalar)
The linear capillary pressure-saturation curve.
Definition: FluidMatrixInteractionBlackoil.hpp:115
static void S(SatContainerT &saturations, const Params &params, const pcContainerT &pc, Scalar)
The saturation-capillary pressure curve.
Definition: FluidMatrixInteractionBlackoil.hpp:139
ParamsT Params
Definition: FluidMatrixInteractionBlackoil.hpp:101
Definition: BlackoilFluid.hpp:32