PolymerProperties.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2012 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_POLYMERPROPERTIES_HEADER_INCLUDED
21 #define OPM_POLYMERPROPERTIES_HEADER_INCLUDED
22 
23 #include <opm/parser/eclipse/Deck/Deck.hpp>
24 #include <opm/parser/eclipse/EclipseState/EclipseState.hpp>
25 
26 #include <cmath>
27 #include <vector>
28 #include <opm/common/ErrorMacros.hpp>
29 
30 
31 namespace Opm
32 {
33 
35  {
36  public:
38  {
39  }
40 
42 
56  PolymerProperties(double c_max,
57  double mix_param,
58  double rock_density,
59  double dead_pore_vol,
60  double res_factor,
61  double c_max_ads,
62  AdsorptionBehaviour ads_index,
63  const std::vector<double>& c_vals_visc,
64  const std::vector<double>& visc_mult_vals,
65  const std::vector<double>& c_vals_ads,
66  const std::vector<double>& ads_vals,
67  const std::vector<double>& water_vel_vals,
68  const std::vector<double>& shear_vrf_vals
69  )
70  : c_max_(c_max),
71  mix_param_(mix_param),
72  rock_density_(rock_density),
73  dead_pore_vol_(dead_pore_vol),
74  res_factor_(res_factor),
75  c_max_ads_(c_max_ads),
76  ads_index_(ads_index),
77  c_vals_visc_(c_vals_visc),
78  visc_mult_vals_(visc_mult_vals),
79  c_vals_ads_(c_vals_ads),
80  ads_vals_(ads_vals),
81  water_vel_vals_(water_vel_vals),
82  shear_vrf_vals_(shear_vrf_vals)
83  {
84  }
85 
86  PolymerProperties(Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclipseState)
87  {
88  readFromDeck(deck, eclipseState);
89  }
90 
91  void set(double c_max,
92  double mix_param,
93  double rock_density,
94  double dead_pore_vol,
95  double res_factor,
96  double c_max_ads,
97  AdsorptionBehaviour ads_index,
98  const std::vector<double>& c_vals_visc,
99  const std::vector<double>& visc_mult_vals,
100  const std::vector<double>& c_vals_ads,
101  const std::vector<double>& ads_vals,
102  const std::vector<double>& water_vel_vals,
103  const std::vector<double>& shear_vrf_vals
104  )
105  {
106  c_max_ = c_max;
107  mix_param_ = mix_param;
108  rock_density_ = rock_density;
109  dead_pore_vol_ = dead_pore_vol;
110  res_factor_ = res_factor;
111  c_max_ads_ = c_max_ads;
112  c_vals_visc_ = c_vals_visc;
113  visc_mult_vals_ = visc_mult_vals;
114  c_vals_ads_ = c_vals_ads;
115  ads_vals_ = ads_vals;
116  ads_index_ = ads_index;
117  water_vel_vals_ = water_vel_vals;
118  shear_vrf_vals_ = shear_vrf_vals;
119  }
120 
121  void readFromDeck(Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclipseState)
122  {
123  // We assume NTMISC=1
124  auto tables = eclipseState->getTableManager();
125  const auto& plymaxTable = tables->getPlymaxTables().getTable<PlymaxTable>(0);
126  const auto plmixparRecord = deck->getKeyword("PLMIXPAR")->getRecord(0);
127 
128  // We also assume that each table has exactly one row...
129  assert(plymaxTable.numRows() == 1);
130 
131  c_max_ = plymaxTable.getPolymerConcentrationColumn()[0];
132  mix_param_ = plmixparRecord->getItem("TODD_LONGSTAFF")->getSIDouble(0);
133 
134  // We assume NTSFUN=1
135  const auto& plyrockTable = tables->getPlyrockTables().getTable<PlyrockTable>(0);
136 
137  // We also assume that each table has exactly one row...
138  assert(plyrockTable.numRows() == 1);
139 
140  dead_pore_vol_ = plyrockTable.getDeadPoreVolumeColumn()[0];
141  res_factor_ = plyrockTable.getResidualResistanceFactorColumn()[0];
142  rock_density_ = plyrockTable.getRockDensityFactorColumn()[0];
143  ads_index_ = static_cast<AdsorptionBehaviour>(plyrockTable.getAdsorbtionIndexColumn()[0]);
144  c_max_ads_ = plyrockTable.getMaxAdsorbtionColumn()[0];
145 
146  // We assume NTPVT=1
147  const auto& plyviscTable = tables->getPlyviscTables().getTable<PlyviscTable>(0);
148 
149 
150  c_vals_visc_ = plyviscTable.getPolymerConcentrationColumn();
151  visc_mult_vals_ = plyviscTable.getViscosityMultiplierColumn();
152 
153  // We assume NTSFUN=1
154  const auto& plyadsTable = tables->getPlyadsTables().getTable<PlyadsTable>(0);
155 
156  c_vals_ads_ = plyadsTable.getPolymerConcentrationColumn();
157  ads_vals_ = plyadsTable.getAdsorbedPolymerColumn();
158 
159  has_plyshlog_ = deck->hasKeyword("PLYSHLOG");
160  has_shrate_ = deck->hasKeyword("SHRATE");
161 
162  if (has_plyshlog_) {
163  // Assuming NTPVT == 1 always
164  const auto& plyshlogTable = tables->getPlyshlogTables().getTable<PlyshlogTable>(0);
165 
166  water_vel_vals_ = plyshlogTable.getWaterVelocityColumn();
167  shear_vrf_vals_ = plyshlogTable.getShearMultiplierColumn();
168 
169  // do the unit version here for the water_vel_vals_
170  Opm::UnitSystem unitSystem = *deck->getActiveUnitSystem();
171  double siFactor;
172  if (has_shrate_) {
173  siFactor = unitSystem.parse("1/Time")->getSIScaling();
174  DeckKeywordConstPtr shrateKeyword = deck->getKeyword("SHRATE");
175  std::vector<double> shrate_readin = shrateKeyword->getSIDoubleData();
176  if (shrate_readin.size() == 1) {
177  shrate_ = shrate_readin[0];
178  } else if (shrate_readin.size() == 0) {
179  shrate_ = 4.8; // default value
180  } else {
181  OPM_THROW(std::logic_error, "Only NTPVT == 1 is allowed for SHRATE keyword now !\n");
182  }
183  } else {
184  siFactor = unitSystem.parse("Length/Time")->getSIScaling();
185  }
186 
187  for (size_t i = 0; i < water_vel_vals_.size(); ++i) {
188  water_vel_vals_[i] *= siFactor;
189  }
190 
191 
192  plyshlog_ref_conc_ = plyshlogTable.getRefPolymerConcentration();
193 
194  if (plyshlogTable.hasRefSalinity()) {
195  has_plyshlog_ref_salinity_ = true;
196  plyshlog_ref_salinity_ = plyshlogTable.getRefSalinity();
197  } else {
198  has_plyshlog_ref_salinity_ = false;
199  }
200 
201  if (plyshlogTable.hasRefTemperature()) {
202  has_plyshlog_ref_temp_ = true;
203  plyshlog_ref_temp_ = plyshlogTable.getRefTemperature();
204  } else {
205  has_plyshlog_ref_temp_ = false;
206  }
207  }
208  }
209 
210  double cMax() const;
211 
212  double mixParam() const;
213 
214  double rockDensity() const;
215 
216  double deadPoreVol() const;
217 
218  double resFactor() const;
219 
220  double cMaxAds() const;
221 
222  int adsIndex() const;
223 
225  bool hasPlyshlog() const;
226 
228  const std::vector<double>& shearWaterVelocity() const;
229 
231  const std::vector<double>& shearViscosityReductionFactor() const;
232 
234  double plyshlogRefConc() const;
235 
237  bool hasPlyshlogRefSalinity() const;
238 
240  bool hasPlyshlogRefTemp() const;
241 
243  double plyshlogRefSalinity() const;
244 
246  double plyshlogRefTemp() const;
247 
249  bool hasShrate() const;
250 
252  double shrate() const;
253 
254  double shearVrf(const double velocity) const;
255 
256  double shearVrfWithDer(const double velocity, double& der) const;
257 
258  double viscMult(double c) const;
259 
260  double viscMultWithDer(double c, double* der) const;
261 
262  void simpleAdsorption(double c, double& c_ads) const;
263 
264  void simpleAdsorptionWithDer(double c, double& c_ads,
265  double& dc_ads_dc) const;
266 
267  void adsorption(double c, double cmax, double& c_ads) const;
268 
269  void adsorptionWithDer(double c, double cmax,
270  double& c_ads, double& dc_ads_dc) const;
271 
272  void effectiveVisc(const double c, const double* visc,
273  double& mu_w_eff) const;
274 
275  void effectiveViscWithDer(const double c, const double* visc
276  , double& mu_w_eff
277  , double dmu_w_eff_dc) const;
278 
279  void effectiveInvVisc(const double c, const double* visc,
280  double& inv_mu_w_eff) const;
281 
282  void effectiveInvViscWithDer(const double c,
283  const double* visc,
284  double& inv_mu_w_eff,
285  double& dinv_mu_w_eff_dc) const;
286  void effectiveRelperm(const double c,
287  const double cmax,
288  const double* relperm,
289  double& eff_relperm_wat) const;
290 
291  void effectiveRelpermWithDer (const double c,
292  const double cmax,
293  const double* relperm,
294  const double* drelperm_ds,
295  double& eff_relperm_wat,
296  double& deff_relperm_wat_ds,
297  double& deff_relperm_wat_dc) const;
298 
299  void effectiveMobilities(const double c,
300  const double cmax,
301  const double* visc,
302  const double* relperm,
303  double* mob) const;
304 
305  void effectiveMobilitiesWithDer(const double c,
306  const double cmax,
307  const double* visc,
308  const double* relperm,
309  const double* drelpermds,
310  double* mob,
311  double* dmob_ds,
312  double& dmobwatdc) const;
313 
314  void effectiveMobilitiesBoth(const double c,
315  const double cmax,
316  const double* visc,
317  const double* relperm,
318  const double* drelperm_ds,
319  double* mob,
320  double* dmob_ds,
321  double& dmobwat_dc,
322  bool if_with_der) const;
323 
324  void effectiveTotalMobility(const double c,
325  const double cmax,
326  const double* visc,
327  const double* relperm,
328  double& totmob) const;
329 
330  void effectiveTotalMobilityWithDer(const double c,
331  const double cmax,
332  const double* visc,
333  const double* relperm,
334  const double* drelpermds,
335  double& totmob,
336  double* dtotmob_dsdc) const;
337 
338  void effectiveTotalMobilityBoth(const double c,
339  const double cmax,
340  const double* visc,
341  const double* relperm,
342  const double* drelperm_ds,
343  double& totmob,
344  double* dtotmob_dsdc,
345  bool if_with_der) const;
346 
347  void computeMc(const double& c, double& mc) const;
348 
349  void computeMcWithDer(const double& c, double& mc,
350  double& dmc_dc) const;
351 
352  void computeMcBoth(const double& c, double& mc,
353  double& dmc_dc, bool if_with_der) const;
354 
356  bool computeShearMultLog(std::vector<double>& water_vel, std::vector<double>& visc_mult, std::vector<double>& shear_mult) const;
357 
358  private:
359  double c_max_;
360  double mix_param_;
361  double rock_density_;
362  double dead_pore_vol_;
363  double res_factor_;
364  double c_max_ads_;
365 
366  bool has_plyshlog_;
367  bool has_shrate_;
368  // Assuming NTPVT == 1 always due to the limitation of the parser
369  // only one SHRATE value
370  // TODO: to be extended later when parser is improved.
371  double shrate_;
372  AdsorptionBehaviour ads_index_;
373  std::vector<double> c_vals_visc_;
374  std::vector<double> visc_mult_vals_;
375  std::vector<double> c_vals_ads_;
376  std::vector<double> ads_vals_;
377  std::vector<double> water_vel_vals_;
378  std::vector<double> shear_vrf_vals_;
379 
380  double plyshlog_ref_conc_;
381  double plyshlog_ref_salinity_;
382  double plyshlog_ref_temp_;
383  bool has_plyshlog_ref_salinity_;
384  bool has_plyshlog_ref_temp_;
385 
386 
387  void simpleAdsorptionBoth(double c, double& c_ads,
388  double& dc_ads_dc, bool if_with_der) const;
389  void adsorptionBoth(double c, double cmax,
390  double& c_ads, double& dc_ads_dc,
391  bool if_with_der) const;
392  void effectiveInvViscBoth(const double c, const double* visc,
393  double& inv_mu_w_eff,
394  double& dinv_mu_w_eff_dc, bool if_with_der) const;
395  void effectiveRelpermBoth(const double c,
396  const double cmax,
397  const double* relperm,
398  const double* drelperm_ds,
399  double& eff_relperm_wat,
400  double& deff_relperm_wat_ds,
401  double& deff_relperm_wat_dc,
402  bool if_with_der) const;
403 
404  };
405 
406 } // namespace Opm
407 
408 #endif // OPM_POLYMERPROPERTIES_HEADER_INCLUDED
void effectiveRelpermWithDer(const double c, const double cmax, const double *relperm, const double *drelperm_ds, double &eff_relperm_wat, double &deff_relperm_wat_ds, double &deff_relperm_wat_dc) const
PolymerProperties(Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclipseState)
Definition: PolymerProperties.hpp:86
double deadPoreVol() const
bool hasPlyshlogRefTemp() const
indicate whether reference temperature is specified in PLYSHLOG
void simpleAdsorptionWithDer(double c, double &c_ads, double &dc_ads_dc) const
void effectiveRelperm(const double c, const double cmax, const double *relperm, double &eff_relperm_wat) const
void computeMcWithDer(const double &c, double &mc, double &dmc_dc) const
Definition: CompressibleTpfaPolymer.hpp:32
void effectiveMobilities(const double c, const double cmax, const double *visc, const double *relperm, double *mob) const
std::vector< double > & cmax
Definition: GravityColumnSolverPolymer_impl.hpp:78
double resFactor() const
double shrate() const
the value of SHRATE
Definition: PolymerProperties.hpp:41
double mixParam() const
void effectiveMobilitiesWithDer(const double c, const double cmax, const double *visc, const double *relperm, const double *drelpermds, double *mob, double *dmob_ds, double &dmobwatdc) const
void effectiveTotalMobilityWithDer(const double c, const double cmax, const double *visc, const double *relperm, const double *drelpermds, double &totmob, double *dtotmob_dsdc) const
double plyshlogRefTemp() const
the reference temperature in PLYSHLOG
bool hasShrate() const
indicate whether SHRATE keyword is specified
double rockDensity() const
void effectiveVisc(const double c, const double *visc, double &mu_w_eff) const
void effectiveViscWithDer(const double c, const double *visc, double &mu_w_eff, double dmu_w_eff_dc) const
Definition: PolymerProperties.hpp:34
PolymerProperties()
Definition: PolymerProperties.hpp:37
double shearVrf(const double velocity) const
void adsorptionWithDer(double c, double cmax, double &c_ads, double &dc_ads_dc) const
void readFromDeck(Opm::DeckConstPtr deck, Opm::EclipseStateConstPtr eclipseState)
Definition: PolymerProperties.hpp:121
bool computeShearMultLog(std::vector< double > &water_vel, std::vector< double > &visc_mult, std::vector< double > &shear_mult) const
Computing the shear multiplier based on the water velocity/shear rate with PLYSHLOG keyword...
double viscMultWithDer(double c, double *der) const
void adsorption(double c, double cmax, double &c_ads) const
void effectiveTotalMobility(const double c, const double cmax, const double *visc, const double *relperm, double &totmob) const
double cMax() const
double shearVrfWithDer(const double velocity, double &der) const
Definition: PolymerProperties.hpp:41
void effectiveTotalMobilityBoth(const double c, const double cmax, const double *visc, const double *relperm, const double *drelperm_ds, double &totmob, double *dtotmob_dsdc, bool if_with_der) const
void effectiveInvVisc(const double c, const double *visc, double &inv_mu_w_eff) const
double plyshlogRefConc() const
the reference polymer concentration in PLYSHLOG
const std::vector< double > & shearViscosityReductionFactor() const
the viscosity reduction factor PLYSHLOG table
const std::vector< double > & shearWaterVelocity() const
the water velocity or water shear rate in PLYSHLOG table
void computeMc(const double &c, double &mc) const
void effectiveMobilitiesBoth(const double c, const double cmax, const double *visc, const double *relperm, const double *drelperm_ds, double *mob, double *dmob_ds, double &dmobwat_dc, bool if_with_der) const
PolymerProperties(double c_max, double mix_param, double rock_density, double dead_pore_vol, double res_factor, double c_max_ads, AdsorptionBehaviour ads_index, const std::vector< double > &c_vals_visc, const std::vector< double > &visc_mult_vals, const std::vector< double > &c_vals_ads, const std::vector< double > &ads_vals, const std::vector< double > &water_vel_vals, const std::vector< double > &shear_vrf_vals)
Definition: PolymerProperties.hpp:56
void simpleAdsorption(double c, double &c_ads) const
void effectiveInvViscWithDer(const double c, const double *visc, double &inv_mu_w_eff, double &dinv_mu_w_eff_dc) const
double cMaxAds() const
double plyshlogRefSalinity() const
the reference salinity in PLYSHLOG
bool hasPlyshlog() const
indicate whehter PLYSHLOG is specified
void computeMcBoth(const double &c, double &mc, double &dmc_dc, bool if_with_der) const
void set(double c_max, double mix_param, double rock_density, double dead_pore_vol, double res_factor, double c_max_ads, AdsorptionBehaviour ads_index, const std::vector< double > &c_vals_visc, const std::vector< double > &visc_mult_vals, const std::vector< double > &c_vals_ads, const std::vector< double > &ads_vals, const std::vector< double > &water_vel_vals, const std::vector< double > &shear_vrf_vals)
Definition: PolymerProperties.hpp:91
double viscMult(double c) const
bool hasPlyshlogRefSalinity() const
indicate wheter reference salinity is specified in PLYSHLOG
AdsorptionBehaviour
Definition: PolymerProperties.hpp:41