opm-common
Loading...
Searching...
No Matches
ConstantCompressibilityBrinePvt.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 2 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*/
27#ifndef OPM_CONSTANT_COMPRESSIBILITY_BRINE_PVT_HPP
28#define OPM_CONSTANT_COMPRESSIBILITY_BRINE_PVT_HPP
29
32
33#include <cstddef>
34#include <vector>
35
36namespace Opm {
37
38template <class Scalar, bool enableThermal, bool enableBrine>
40
41#if HAVE_ECL_INPUT
42class EclipseState;
43class Schedule;
44#endif
45
50template <class Scalar>
52{
53public:
54 using TabulatedFunction = Tabulated1DFunction<Scalar>;
55
56#if HAVE_ECL_INPUT
61 void initFromState(const EclipseState& eclState, const Schedule&);
62#endif
63
64 void setNumRegions(std::size_t numRegions);
65
66 void setVapPars(const Scalar, const Scalar)
67 {
68 }
69
73 void setReferenceDensities(unsigned regionIdx,
74 Scalar /*rhoRefOil*/,
75 Scalar /*rhoRefGas*/,
76 Scalar rhoRefWater)
77 { waterReferenceDensity_[regionIdx] = rhoRefWater; }
78
82 void initEnd()
83 { }
84
88 unsigned numRegions() const
89 { return waterReferenceDensity_.size(); }
90
94 template <class Evaluation>
95 Evaluation internalEnergy(unsigned,
96 const Evaluation&,
97 const Evaluation&,
98 const Evaluation&,
99 const Evaluation&) const
100 {
101 throw std::runtime_error("Requested the enthalpy of water but the thermal option is not enabled");
102 }
103
104 Scalar hVap(unsigned) const
105 {
106 throw std::runtime_error("Requested the hvap of oil but the thermal option is not enabled");
107 }
108
112 template <class Evaluation>
113 Evaluation viscosity(unsigned regionIdx,
114 const Evaluation& temperature,
115 const Evaluation& pressure,
116 const Evaluation& Rsw,
117 const Evaluation& saltconcentration) const
118 {
119 // cf. ECLiPSE 2013.2 technical description, p. 114
120 Scalar pRef = referencePressure_[regionIdx];
121 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
122 const Evaluation Cv = viscosibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
123 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
124 const Evaluation Y = (C-Cv)* (pressure - pRef);
125 Evaluation MuwRef = viscosityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
126
127 const Evaluation& bw = inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rsw, saltconcentration);
128
129 return MuwRef * BwRef * bw / (1 + Y * (1 + Y/2));
130 }
131
132
136 template <class Evaluation>
137 Evaluation saturatedViscosity(unsigned regionIdx,
138 const Evaluation& temperature,
139 const Evaluation& pressure,
140 const Evaluation& saltconcentration) const
141 {
142 Scalar pRef = referencePressure_[regionIdx];
143 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
144 const Evaluation Cv = viscosibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
145 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
146 const Evaluation Y = (C-Cv)* (pressure - pRef);
147 Evaluation MuwRef = viscosityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
148
149 const Evaluation& bw = saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure, saltconcentration);
150
151 return MuwRef * BwRef * bw / (1 + Y * (1 + Y/2));
152 }
153
157 template <class Evaluation>
158 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
159 const Evaluation& temperature,
160 const Evaluation& pressure,
161 const Evaluation& saltconcentration) const
162 {
163 Evaluation Rsw = 0.0;
164 return inverseFormationVolumeFactor(regionIdx, temperature, pressure,
165 Rsw, saltconcentration);
166 }
167
170 template <class Evaluation>
171 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
172 const Evaluation& /*temperature*/,
173 const Evaluation& pressure,
174 const Evaluation& /*Rsw*/,
175 const Evaluation& saltconcentration) const
176 {
177 Scalar pRef = referencePressure_[regionIdx];
178
179 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
180 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
181 const Evaluation X = C * (pressure - pRef);
182
183 return (1.0 + X * (1.0 + X / 2.0)) / BwRef;
184
185 }
186
190 template <class FluidState, class LhsEval = typename FluidState::Scalar>
191 std::pair<LhsEval, LhsEval>
192 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
193 {
194 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(FluidState::waterPhaseIdx));
195 const LhsEval& saltConcentration
196 = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
197 const auto segIdx = this->formationVolumeTables_[regionIdx]
198 .findSegmentIndex(saltConcentration, /*extrapolate=*/ true);
199
200 // Calculate bw.
201 const Scalar pRef = referencePressure_[regionIdx];
202 const LhsEval BwRef = formationVolumeTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
203 const LhsEval C = compressibilityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
204 const LhsEval X = C * (pressure - pRef);
205 const LhsEval bw = (1.0 + X * (1.0 + X / 2.0)) / BwRef;
206
207 // Calculate muw.
208 const LhsEval Cv = viscosibilityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
209 const LhsEval Y = (C - Cv) * (pressure - pRef);
210 const LhsEval MuwRef = viscosityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
211 const LhsEval muw = MuwRef * BwRef * bw / (1.0 + Y * (1.0 + Y / 2.0));
212
213 return { bw, muw };
214 }
215
222 template <class Evaluation>
223 Evaluation saturationPressure(unsigned /*regionIdx*/,
224 const Evaluation& /*temperature*/,
225 const Evaluation& /*Rs*/,
226 const Evaluation& /*saltconcentration*/) const
227 { return 0.0; /* this is dead water, so there isn't any meaningful saturation pressure! */ }
228
232 template <class Evaluation>
233 Evaluation saturatedGasDissolutionFactor(unsigned /*regionIdx*/,
234 const Evaluation& /*temperature*/,
235 const Evaluation& /*pressure*/,
236 const Evaluation& /*saltconcentration*/) const
237 { return 0.0; /* this is dead water! */ }
238
239 template <class Evaluation>
240 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
241 const Evaluation& /*pressure*/,
242 unsigned /*compIdx*/) const
243 {
244 throw std::runtime_error("Not implemented: The PVT model does not provide "
245 "a diffusionCoefficient()");
246 }
247
248 Scalar waterReferenceDensity(unsigned regionIdx) const
249 { return waterReferenceDensity_[regionIdx]; }
250
251 const std::vector<Scalar>& referencePressure() const
252 { return referencePressure_; }
253
254 const std::vector<TabulatedFunction>& formationVolumeTables() const
255 { return formationVolumeTables_; }
256
257 const std::vector<TabulatedFunction>& compressibilityTables() const
258 { return compressibilityTables_; }
259
260 const std::vector<TabulatedFunction>& viscosityTables() const
261 { return viscosityTables_; }
262
263 const std::vector<TabulatedFunction>& viscosibilityTables() const
264 { return viscosibilityTables_; }
265
266private:
267 std::vector<Scalar> waterReferenceDensity_{};
268 std::vector<Scalar> referencePressure_{};
269 std::vector<TabulatedFunction> formationVolumeTables_{};
270 std::vector<TabulatedFunction> compressibilityTables_{};
271 std::vector<TabulatedFunction> viscosityTables_{};
272 std::vector<TabulatedFunction> viscosibilityTables_{};
273};
274
275} // namespace Opm
276
277#endif
Implements a linearly interpolated scalar function that depends on one variable.
This class represents the Pressure-Volume-Temperature relations of the gas phase without vaporized oi...
Definition ConstantCompressibilityBrinePvt.hpp:52
Evaluation saturatedViscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition ConstantCompressibilityBrinePvt.hpp:137
Evaluation viscosity(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &Rsw, const Evaluation &saltconcentration) const
Returns the dynamic viscosity [Pa s] of the fluid phase given a set of parameters.
Definition ConstantCompressibilityBrinePvt.hpp:113
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition ConstantCompressibilityBrinePvt.hpp:88
Evaluation inverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &, const Evaluation &pressure, const Evaluation &, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition ConstantCompressibilityBrinePvt.hpp:171
void setReferenceDensities(unsigned regionIdx, Scalar, Scalar, Scalar rhoRefWater)
Set the water reference density [kg / m^3].
Definition ConstantCompressibilityBrinePvt.hpp:73
Evaluation saturationPressure(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the saturation pressure of the water phase [Pa] depending on its mass fraction of the gas com...
Definition ConstantCompressibilityBrinePvt.hpp:223
std::pair< LhsEval, LhsEval > inverseFormationVolumeFactorAndViscosity(const FluidState &fluidState, unsigned regionIdx)
Returns the formation volume factor [-] and viscosity [Pa s] of the fluid phase.
Definition ConstantCompressibilityBrinePvt.hpp:192
Evaluation saturatedGasDissolutionFactor(unsigned, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the gas dissolution factor [m^3/m^3] of the water phase.
Definition ConstantCompressibilityBrinePvt.hpp:233
Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx, const Evaluation &temperature, const Evaluation &pressure, const Evaluation &saltconcentration) const
Returns the formation volume factor [-] of the fluid phase.
Definition ConstantCompressibilityBrinePvt.hpp:158
void initEnd()
Finish initializing the water phase PVT properties.
Definition ConstantCompressibilityBrinePvt.hpp:82
Evaluation internalEnergy(unsigned, const Evaluation &, const Evaluation &, const Evaluation &, const Evaluation &) const
Returns the specific enthalpy [J/kg] of water given a set of parameters.
Definition ConstantCompressibilityBrinePvt.hpp:95
Definition EclipseState.hpp:62
Definition Schedule.hpp:101
Implements a linearly interpolated scalar function that depends on one variable.
Definition Tabulated1DFunction.hpp:51
This class represents the Pressure-Volume-Temperature relations of the water phase in the black-oil m...
Definition WaterPvtMultiplexer.hpp:90
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Tabulated1DFunction.hpp:41