20 #ifndef OPM_PVTCONSTCOMPR_HEADER_INCLUDED
21 #define OPM_PVTCONSTCOMPR_HEADER_INCLUDED
24 #include <opm/common/ErrorMacros.hpp>
52 int numRegions = pvtwKeyword->size();
54 ref_press_.resize(numRegions);
55 ref_B_.resize(numRegions);
56 comp_.resize(numRegions);
57 viscosity_.resize(numRegions);
58 visc_comp_.resize(numRegions);
60 for (
int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
61 Opm::DeckRecordConstPtr pvtwRecord = pvtwKeyword->getRecord(regionIdx);
63 ref_press_[regionIdx] = pvtwRecord->getItem(
"P_REF")->getSIDouble(0);
64 ref_B_[regionIdx] = pvtwRecord->getItem(
"WATER_VOL_FACTOR")->getSIDouble(0);
65 comp_[regionIdx] = pvtwRecord->getItem(
"WATER_COMPRESSIBILITY")->getSIDouble(0);
66 viscosity_[regionIdx] = pvtwRecord->getItem(
"WATER_VISCOSITY")->getSIDouble(0);
67 visc_comp_[regionIdx] = pvtwRecord->getItem(
"WATER_VISCOSIBILITY")->getSIDouble(0);
73 int numRegions = pvcdoKeyword->size();
75 ref_press_.resize(numRegions);
76 ref_B_.resize(numRegions);
77 comp_.resize(numRegions);
78 viscosity_.resize(numRegions);
79 visc_comp_.resize(numRegions);
81 for (
int regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
82 Opm::DeckRecordConstPtr pvcdoRecord = pvcdoKeyword->getRecord(regionIdx);
84 ref_press_[regionIdx] = pvcdoRecord->getItem(
"P_REF")->getSIDouble(0);
85 ref_B_[regionIdx] = pvcdoRecord->getItem(
"OIL_VOL_FACTOR")->getSIDouble(0);
86 comp_[regionIdx] = pvcdoRecord->getItem(
"OIL_COMPRESSIBILITY")->getSIDouble(0);
87 viscosity_[regionIdx] = pvcdoRecord->getItem(
"OIL_VISCOSITY")->getSIDouble(0);
88 visc_comp_[regionIdx] = pvcdoRecord->getItem(
"OIL_VISCOSIBILITY")->getSIDouble(0);
109 virtual void mu(
const int n,
110 const int* pvtRegionIdx,
114 double* output_mu)
const
117 for (
int i = 0; i < n; ++i) {
119 int tableIdx = getTableIndex_(pvtRegionIdx, i);
120 double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
121 output_mu[i] = viscosity_[tableIdx]/(1.0 + x + 0.5*x*x);
125 virtual void mu(
const int n,
126 const int* pvtRegionIdx,
131 double* output_dmudp,
132 double* output_dmudr)
const
135 for (
int i = 0; i < n; ++i) {
137 int tableIdx = getTableIndex_(pvtRegionIdx, i);
138 double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
139 double d = (1.0 + x + 0.5*x*x);
140 output_mu[i] = viscosity_[tableIdx]/d;
141 output_dmudp[i] = (viscosity_[tableIdx]/(d*d))*(1+x) * visc_comp_[tableIdx];
143 std::fill(output_dmudr, output_dmudr + n, 0.0);
146 virtual void mu(
const int n,
147 const int* pvtRegionIdx,
153 double* output_dmudp,
154 double* output_dmudr)
const
157 for (
int i = 0; i < n; ++i) {
159 int tableIdx = getTableIndex_(pvtRegionIdx, i);
160 double x = -visc_comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
161 double d = (1.0 + x + 0.5*x*x);
162 output_mu[i] = viscosity_[tableIdx]/d;
163 output_dmudp[i] = (viscosity_[tableIdx]/(d*d))*(1+x) * visc_comp_[tableIdx];
165 std::fill(output_dmudr, output_dmudr + n, 0.0);
168 virtual void B(
const int n,
169 const int* pvtRegionIdx,
173 double* output_B)
const
176 for (
int i = 0; i < n; ++i) {
178 int tableIdx = getTableIndex_(pvtRegionIdx, i);
179 double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
180 output_B[i] = ref_B_[tableIdx]/(1.0 + x + 0.5*x*x);
185 const int* pvtRegionIdx,
190 double* output_dBdp)
const
193 for (
int i = 0; i < n; ++i) {
194 int tableIdx = getTableIndex_(pvtRegionIdx, i);
195 double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
196 double d = (1.0 + x + 0.5*x*x);
197 output_B[i] = ref_B_[tableIdx]/d;
198 output_dBdp[i] = (-ref_B_[tableIdx]/(d*d))*(1 + x) * comp_[tableIdx];
202 virtual void b(
const int n,
203 const int* pvtRegionIdx,
209 double* output_dbdr)
const
212 for (
int i = 0; i < n; ++i) {
214 int tableIdx = getTableIndex_(pvtRegionIdx, i);
215 double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
216 double d = (1.0 + x + 0.5*x*x);
219 output_b[i] = d/ref_B_[tableIdx];
220 output_dbdp[i] = (1 + x) * comp_[tableIdx]/ref_B_[tableIdx];
223 std::fill(output_dbdr, output_dbdr + n, 0.0);
226 virtual void b(
const int n,
227 const int* pvtRegionIdx,
234 double* output_dbdr)
const
237 for (
int i = 0; i < n; ++i) {
239 int tableIdx = getTableIndex_(pvtRegionIdx, i);
240 double x = comp_[tableIdx]*(p[i] - ref_press_[tableIdx]);
241 double d = (1.0 + x + 0.5*x*x);
244 output_b[i] = d/ref_B_[tableIdx];
245 output_dbdp[i] = (1 + x) * comp_[tableIdx]/ref_B_[tableIdx];
247 std::fill(output_dbdr, output_dbdr + n, 0.0);
253 double* output_rsSat,
254 double* output_drsSatdp)
const
256 std::fill(output_rsSat, output_rsSat + n, 0.0);
257 std::fill(output_drsSatdp, output_drsSatdp + n, 0.0);
263 double* output_rvSat,
264 double* output_drvSatdp)
const
266 std::fill(output_rvSat, output_rvSat + n, 0.0);
267 std::fill(output_drvSatdp, output_drvSatdp + n, 0.0);
270 virtual void R(
const int n,
274 double* output_R)
const
276 std::fill(output_R, output_R + n, 0.0);
284 double* output_dRdp)
const
286 std::fill(output_R, output_R + n, 0.0);
287 std::fill(output_dRdp, output_dRdp + n, 0.0);
291 int getTableIndex_(
const int* pvtTableIdx,
int cellIdx)
const
295 return pvtTableIdx[cellIdx];
300 std::vector<double> ref_press_;
301 std::vector<double> ref_B_;
302 std::vector<double> comp_;
303 std::vector<double> viscosity_;
304 std::vector<double> visc_comp_;
310 #endif // OPM_PVTCONSTCOMPR_HEADER_INCLUDED
Definition: BlackoilPhases.hpp:48
Definition: AnisotropicEikonal.hpp:43
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, const PhasePresence *, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: PvtConstCompr.hpp:146
virtual void dRdp(const int n, const int *, const double *, const double *, double *output_R, double *output_dRdp) const
Solution factor and p-derivative as functions of p and z.
Definition: PvtConstCompr.hpp:279
virtual void R(const int n, const int *, const double *, const double *, double *output_R) const
Solution factor as a function of p and z.
Definition: PvtConstCompr.hpp:270
PvtConstCompr()
Definition: PvtConstCompr.hpp:47
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: PvtConstCompr.hpp:202
virtual void dBdp(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_B, double *output_dBdp) const
Formation volume factor and p-derivative as functions of p and z.
Definition: PvtConstCompr.hpp:184
void initFromOil(Opm::DeckKeywordConstPtr pvcdoKeyword)
Definition: PvtConstCompr.hpp:71
virtual void B(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_B) const
Formation volume factor as a function of p, T and z.
Definition: PvtConstCompr.hpp:168
virtual void b(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, const PhasePresence *, double *output_b, double *output_dbdp, double *output_dbdr) const
Definition: PvtConstCompr.hpp:226
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_mu) const
Viscosity as a function of p, T and z.
Definition: PvtConstCompr.hpp:109
Definition: PvtConstCompr.hpp:44
void initFromWater(Opm::DeckKeywordConstPtr pvtwKeyword)
Definition: PvtConstCompr.hpp:50
Definition: PvtInterface.hpp:30
virtual void rsSat(const int n, const int *, const double *, double *output_rsSat, double *output_drsSatdp) const
Solution gas/oil ratio and its derivatives at saturated conditions as a function of p...
Definition: PvtConstCompr.hpp:250
virtual void rvSat(const int n, const int *, const double *, double *output_rvSat, double *output_drvSatdp) const
Vapor oil/gas ratio and its derivatives at saturated conditions as a function of p.
Definition: PvtConstCompr.hpp:260
PvtConstCompr(double visc)
Create a PVT object with a given viscosity that assumes all fluid phases to be incompressible.
Definition: PvtConstCompr.hpp:96
const double visc
Definition: Units.hpp:152
virtual void mu(const int n, const int *pvtRegionIdx, const double *p, const double *, const double *, double *output_mu, double *output_dmudp, double *output_dmudr) const
Definition: PvtConstCompr.hpp:125
virtual ~PvtConstCompr()
Definition: PvtConstCompr.hpp:105