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
41class EclipseState;
42class Schedule;
43
48template <class Scalar>
50{
51public:
52 using TabulatedFunction = Tabulated1DFunction<Scalar>;
53
58 void initFromState(const EclipseState& eclState, const Schedule&);
59
60 void setNumRegions(std::size_t numRegions);
61
62 void setVapPars(const Scalar, const Scalar)
63 {
64 }
65
69 void setReferenceDensities(unsigned regionIdx,
70 Scalar /*rhoRefOil*/,
71 Scalar /*rhoRefGas*/,
72 Scalar rhoRefWater)
73 { waterReferenceDensity_[regionIdx] = rhoRefWater; }
74
78 void initEnd()
79 { }
80
84 unsigned numRegions() const
85 { return waterReferenceDensity_.size(); }
86
90 template <class Evaluation>
91 Evaluation internalEnergy(unsigned,
92 const Evaluation&,
93 const Evaluation&,
94 const Evaluation&,
95 const Evaluation&) const
96 {
97 throw std::runtime_error("Requested the enthalpy of water but the thermal option is not enabled");
98 }
99
100 Scalar hVap(unsigned) const
101 {
102 throw std::runtime_error("Requested the hvap of oil but the thermal option is not enabled");
103 }
104
108 template <class Evaluation>
109 Evaluation viscosity(unsigned regionIdx,
110 const Evaluation& temperature,
111 const Evaluation& pressure,
112 const Evaluation& Rsw,
113 const Evaluation& saltconcentration) const
114 {
115 // cf. ECLiPSE 2013.2 technical description, p. 114
116 Scalar pRef = referencePressure_[regionIdx];
117 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
118 const Evaluation Cv = viscosibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
119 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
120 const Evaluation Y = (C-Cv)* (pressure - pRef);
121 Evaluation MuwRef = viscosityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
122
123 const Evaluation& bw = inverseFormationVolumeFactor(regionIdx, temperature, pressure, Rsw, saltconcentration);
124
125 return MuwRef * BwRef * bw / (1 + Y * (1 + Y/2));
126 }
127
128
132 template <class Evaluation>
133 Evaluation saturatedViscosity(unsigned regionIdx,
134 const Evaluation& temperature,
135 const Evaluation& pressure,
136 const Evaluation& saltconcentration) const
137 {
138 Scalar pRef = referencePressure_[regionIdx];
139 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
140 const Evaluation Cv = viscosibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
141 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
142 const Evaluation Y = (C-Cv)* (pressure - pRef);
143 Evaluation MuwRef = viscosityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
144
145 const Evaluation& bw = saturatedInverseFormationVolumeFactor(regionIdx, temperature, pressure, saltconcentration);
146
147 return MuwRef * BwRef * bw / (1 + Y * (1 + Y/2));
148 }
149
153 template <class Evaluation>
154 Evaluation saturatedInverseFormationVolumeFactor(unsigned regionIdx,
155 const Evaluation& temperature,
156 const Evaluation& pressure,
157 const Evaluation& saltconcentration) const
158 {
159 Evaluation Rsw = 0.0;
160 return inverseFormationVolumeFactor(regionIdx, temperature, pressure,
161 Rsw, saltconcentration);
162 }
163
166 template <class Evaluation>
167 Evaluation inverseFormationVolumeFactor(unsigned regionIdx,
168 const Evaluation& /*temperature*/,
169 const Evaluation& pressure,
170 const Evaluation& /*Rsw*/,
171 const Evaluation& saltconcentration) const
172 {
173 Scalar pRef = referencePressure_[regionIdx];
174
175 const Evaluation BwRef = formationVolumeTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
176 const Evaluation C = compressibilityTables_[regionIdx].eval(saltconcentration, /*extrapolate=*/true);
177 const Evaluation X = C * (pressure - pRef);
178
179 return (1.0 + X * (1.0 + X / 2.0)) / BwRef;
180
181 }
182
186 template <class FluidState, class LhsEval = typename FluidState::ValueType>
187 std::pair<LhsEval, LhsEval>
188 inverseFormationVolumeFactorAndViscosity(const FluidState& fluidState, unsigned regionIdx)
189 {
190 const LhsEval& pressure = decay<LhsEval>(fluidState.pressure(FluidState::waterPhaseIdx));
191 const LhsEval& saltConcentration
192 = BlackOil::template getSaltConcentration_<FluidState, LhsEval>(fluidState, regionIdx);
193 const auto segIdx = this->formationVolumeTables_[regionIdx]
194 .findSegmentIndex(saltConcentration, /*extrapolate=*/ true);
195
196 // Calculate bw.
197 const Scalar pRef = referencePressure_[regionIdx];
198 const LhsEval BwRef = formationVolumeTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
199 const LhsEval C = compressibilityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
200 const LhsEval X = C * (pressure - pRef);
201 const LhsEval bw = (1.0 + X * (1.0 + X / 2.0)) / BwRef;
202
203 // Calculate muw.
204 const LhsEval Cv = viscosibilityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
205 const LhsEval Y = (C - Cv) * (pressure - pRef);
206 const LhsEval MuwRef = viscosityTables_[regionIdx].eval(saltConcentration, SegmentIndex{segIdx});
207 const LhsEval muw = MuwRef * BwRef * bw / (1.0 + Y * (1.0 + Y / 2.0));
208
209 return { bw, muw };
210 }
211
218 template <class Evaluation>
219 Evaluation saturationPressure(unsigned /*regionIdx*/,
220 const Evaluation& /*temperature*/,
221 const Evaluation& /*Rs*/,
222 const Evaluation& /*saltconcentration*/) const
223 { return 0.0; /* this is dead water, so there isn't any meaningful saturation pressure! */ }
224
228 template <class Evaluation>
229 Evaluation saturatedGasDissolutionFactor(unsigned /*regionIdx*/,
230 const Evaluation& /*temperature*/,
231 const Evaluation& /*pressure*/,
232 const Evaluation& /*saltconcentration*/) const
233 { return 0.0; /* this is dead water! */ }
234
235 template <class Evaluation>
236 Evaluation diffusionCoefficient(const Evaluation& /*temperature*/,
237 const Evaluation& /*pressure*/,
238 unsigned /*compIdx*/,
239 unsigned /*regionIdx*/ = 0) const
240 {
241 throw std::runtime_error("Not implemented: The PVT model does not provide "
242 "a diffusionCoefficient()");
243 }
244
245 Scalar waterReferenceDensity(unsigned regionIdx) const
246 { return waterReferenceDensity_[regionIdx]; }
247
248 const std::vector<Scalar>& referencePressure() const
249 { return referencePressure_; }
250
251 const std::vector<TabulatedFunction>& formationVolumeTables() const
252 { return formationVolumeTables_; }
253
254 const std::vector<TabulatedFunction>& compressibilityTables() const
255 { return compressibilityTables_; }
256
257 const std::vector<TabulatedFunction>& viscosityTables() const
258 { return viscosityTables_; }
259
260 const std::vector<TabulatedFunction>& viscosibilityTables() const
261 { return viscosibilityTables_; }
262
263private:
264 std::vector<Scalar> waterReferenceDensity_{};
265 std::vector<Scalar> referencePressure_{};
266 std::vector<TabulatedFunction> formationVolumeTables_{};
267 std::vector<TabulatedFunction> compressibilityTables_{};
268 std::vector<TabulatedFunction> viscosityTables_{};
269 std::vector<TabulatedFunction> viscosibilityTables_{};
270};
271
272} // namespace Opm
273
274#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:50
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:133
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:109
unsigned numRegions() const
Return the number of PVT regions which are considered by this PVT-object.
Definition ConstantCompressibilityBrinePvt.hpp:84
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:167
void setReferenceDensities(unsigned regionIdx, Scalar, Scalar, Scalar rhoRefWater)
Set the water reference density [kg / m^3].
Definition ConstantCompressibilityBrinePvt.hpp:69
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:219
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:188
void initFromState(const EclipseState &eclState, const Schedule &)
Sets the pressure-dependent water viscosity and density using a table stemming from the Eclipse PVTWS...
Definition ConstantCompressibilityBrinePvt.cpp:38
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:229
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:154
void initEnd()
Finish initializing the water phase PVT properties.
Definition ConstantCompressibilityBrinePvt.hpp:78
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:91
Definition EclipseState.hpp:66
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:88
This class implements a small container which holds the transmissibility mulitpliers for all the face...
Definition Exceptions.hpp:30
Definition Tabulated1DFunction.hpp:41