RockAnisotropicRelperm.hpp
Go to the documentation of this file.
1 //===========================================================================
2 //
3 // File: RockAnisotropicRelperm.hpp
4 //
5 // Created: Fri Oct 23 14:11:24 2009
6 //
7 // Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8 // Bård Skaflestad <bard.skaflestad@sintef.no>
9 //
10 // $Date$
11 //
12 // $Revision$
13 //
14 //===========================================================================
15 
16 /*
17  Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
18  Copyright 2009, 2010 Statoil ASA.
19 
20  This file is part of The Open Reservoir Simulator Project (OpenRS).
21 
22  OpenRS is free software: you can redistribute it and/or modify
23  it under the terms of the GNU General Public License as published by
24  the Free Software Foundation, either version 3 of the License, or
25  (at your option) any later version.
26 
27  OpenRS is distributed in the hope that it will be useful,
28  but WITHOUT ANY WARRANTY; without even the implied warranty of
29  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  GNU General Public License for more details.
31 
32  You should have received a copy of the GNU General Public License
33  along with OpenRS. If not, see <http://www.gnu.org/licenses/>.
34 */
35 
36 #ifndef OPENRS_ROCKANISOTROPICRELPERM_HEADER
37 #define OPENRS_ROCKANISOTROPICRELPERM_HEADER
38 
39 
40 #include <dune/common/fvector.hh>
41 #include <opm/core/utility/NonuniformTableLinear.hpp>
43 
44 #include <fstream>
45 #include <vector>
46 #include <algorithm>
47 #include <limits>
48 #include <iostream>
49 
50 namespace Opm
51 {
52 
54  {
55  public:
56 
57  void setUseJfunctionScaling(const bool use_j)
58  {
59  if (use_j) {
60  OPM_THROW(std::runtime_error, "RockAnisotropicRelperm cannot use J-scaling.");
61  }
62  }
63 
64  void setSigmaAndTheta(const double, const double)
65  {
66  OPM_THROW(std::runtime_error, "RockAnisotropicRelperm cannot accept sigma and theta arguments.");
67  }
68 
69 
70  template <template <class> class SP, class OP>
71  void kr(const int phase_index, const double saturation, FullMatrix<double, SP, OP>& kr_value) const
72  {
73  zero(kr_value);
74  kr_value(0,0) = krxx_[phase_index](saturation);
75  kr_value(1,1) = kryy_[phase_index](saturation);
76  kr_value(2,2) = krzz_[phase_index](saturation);
77  }
78 
79  template <template <class> class SP, class OP>
80  double capPress(const FullMatrix<double, SP, OP>& /*perm*/, const double /*poro*/, const double saturation) const
81  {
82  return cap_press_(saturation);
83  }
84 
85  template <template <class> class SP, class OP>
86  double satFromCapPress(const FullMatrix<double, SP, OP>& /*perm*/, const double /*poro*/, const double cp) const
87  {
88  return cap_press_.inverse(cp);
89  }
90 
91  void read(const std::string& directory, const std::string& specification)
92  {
93  // For this type of rock, the specification is a line with two file names.
94  std::istringstream specstream(specification);
95  std::string rockname[2];
96  specstream >> rockname[0] >> rockname[1];
97  for (int i = 0; i < 2; ++i) {
98  std::string rockfilename = directory + rockname[i];
99  std::ifstream rock_stream(rockfilename.c_str());
100  if (!rock_stream) {
101  OPM_THROW(std::runtime_error, "Could not open file " << rockfilename);
102  }
103  readAnisoFormat(rock_stream, i);
104  }
105  }
106 
107  private:
108  void readAnisoFormat(std::istream& is, int phase)
109  {
110  // Ignore comments.
111  while (is.peek() == '#') {
112  is.ignore(std::numeric_limits<int>::max(), '\n');
113  }
114  typedef Dune::FieldVector<double, 5> Data;
115  std::istream_iterator<Data> start(is);
116  std::istream_iterator<Data> end;
117  std::vector<Data> data(start, end);
118  if (!is.eof()) {
119  OPM_THROW(std::runtime_error, "Reading stopped but we're not at eof: something went wrong reading data.");
120  }
121  // Making sure that the tables start/end at zero relperm.
122  bool all_ok = true;
123  for (int i = 2; i < 5; ++i) {
124  if (phase == 0 && data[0][i] != 0.0) { // Water table
125  data[0][i] = 0.0;
126  all_ok = false;
127  } else if (phase == 1 && data.back()[i] != 0.0) { // Oil table
128  data.back()[i] = 0.0;
129  all_ok = false;
130  }
131  }
132  if (!all_ok) {
133  std::cout << "Warning: Relperm tables were modified to go to zero at the beginning/end." << std::endl;
134  }
135  std::vector<double> pcow, svals, krxx, kryy, krzz;
136  for (int i = 0; i < int(data.size()); ++i) {
137  pcow.push_back(data[i][0]);
138  svals.push_back(data[i][1]);
139  krxx.push_back(data[i][2]);
140  kryy.push_back(data[i][3]);
141  krzz.push_back(data[i][4]);
142  }
143  krxx_[phase] = TabFunc(svals, krxx);
144  kryy_[phase] = TabFunc(svals, kryy);
145  krzz_[phase] = TabFunc(svals, krzz);
146  if (phase == 0) {
147  cap_press_ = TabFunc(svals, pcow);
148  } else {
149  assert(phase == 1);
150  assert(cap_press_ == TabFunc(svals, pcow));
151  }
152  }
153 
154  typedef NonuniformTableLinear<double> TabFunc;
155  TabFunc krxx_[2];
156  TabFunc kryy_[2];
157  TabFunc krzz_[2];
158  TabFunc cap_press_;
159  };
160 
161 } // namespace Opm
162 
163 
164 
165 #endif // OPENRS_ROCKANISOTROPICRELPERM_HEADER
void setUseJfunctionScaling(const bool use_j)
Definition: RockAnisotropicRelperm.hpp:57
Definition: BlackoilFluid.hpp:31
void read(const std::string &directory, const std::string &specification)
Definition: RockAnisotropicRelperm.hpp:91
Definition: RockAnisotropicRelperm.hpp:53
double satFromCapPress(const FullMatrix< double, SP, OP > &, const double, const double cp) const
Definition: RockAnisotropicRelperm.hpp:86
double capPress(const FullMatrix< double, SP, OP > &, const double, const double saturation) const
Definition: RockAnisotropicRelperm.hpp:80
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603
void kr(const int phase_index, const double saturation, FullMatrix< double, SP, OP > &kr_value) const
Definition: RockAnisotropicRelperm.hpp:71
void setSigmaAndTheta(const double, const double)
Definition: RockAnisotropicRelperm.hpp:64