Rock_impl.hpp
Go to the documentation of this file.
1/*
2 Copyright 2010 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_ROCK_IMPL_HEADER_INCLUDED
21#define OPM_ROCK_IMPL_HEADER_INCLUDED
22
23
24#include <fstream>
25#include <array>
26#include <opm/core/io/eclipse/EclipseGridInspector.hpp>
27
28namespace Opm
29{
30
31 // ----- Methods of Rock -----
32
33
34
35
36 template <int dim>
38 : permeability_kind_(Invalid)
39 {
40 }
41
42
43 template <int dim>
44 void Rock<dim>::init(Opm::DeckConstPtr deck,
45 const std::vector<int>& global_cell,
46 const double perm_threshold)
47 {
48 // This code is mostly copied from ReservoirPropertyCommon::init(...).
49 static_assert(dim == 3, "");
50
51 permfield_valid_.assign(global_cell.size(), false);
52
53 assignPorosity (deck, global_cell);
54 assignPermeability(deck, global_cell, perm_threshold);
55 }
56
57
58 template <int dim>
59 void Rock<dim>::init(const int num_cells,
60 const double uniform_poro,
61 const double uniform_perm)
62 {
63 permfield_valid_.assign(num_cells, true);
64 porosity_.assign(num_cells, uniform_poro);
65 permeability_.assign(dim*dim*num_cells, 0.0);
66 for (int i = 0; i < num_cells; ++i) {
67 SharedPermTensor K = permeabilityModifiable(i);
68 for (int dd = 0; dd < dim; ++dd) {
69 K(dd, dd) = uniform_perm;
70 }
71 }
72 }
73
74
75 template <int dim>
76 double Rock<dim>::porosity(int cell_index) const
77 {
78 return porosity_[cell_index];
79 }
80
81
82 template <int dim>
84 Rock<dim>::permeability(int cell_index) const
85 {
86 assert (permfield_valid_[cell_index]);
87
88 const PermTensor K(dim, dim, &permeability_[dim*dim*cell_index]);
89 return K;
90 }
91
92
93 template <int dim>
96 {
97 // Typically only used for assigning synthetic perm values.
98 SharedPermTensor K(dim, dim, &permeability_[dim*dim*cell_index]);
99
100 // Trust caller!
101 permfield_valid_[cell_index] = std::vector<unsigned char>::value_type(1);
102
103 return K;
104 }
105
106
107
108
109 // ------ Private methods ------
110
111
112
113
114 template <int dim>
115 void Rock<dim>::assignPorosity(Opm::DeckConstPtr deck,
116 const std::vector<int>& global_cell)
117 {
118 porosity_.assign(global_cell.size(), 1.0);
119
120 if (deck->hasKeyword("PORO")) {
121 const std::vector<double>& poro = deck->getKeyword("PORO").getSIDoubleData();
122
123 for (int c = 0; c < int(porosity_.size()); ++c) {
124 porosity_[c] = poro[global_cell[c]];
125 }
126 }
127 }
128
129
130
131 template <int dim>
132 void Rock<dim>::assignPermeability(Opm::DeckConstPtr deck,
133 const std::vector<int>& global_cell,
134 double perm_threshold)
135 {
136 Opm::EclipseGridInspector insp(deck);
137 std::array<int, 3> dims = insp.gridSize();
138 int num_global_cells = dims[0]*dims[1]*dims[2];
139 assert (num_global_cells > 0);
140
141 permeability_.assign(dim * dim * global_cell.size(), 0.0);
142
143 std::vector<const std::vector<double>*> tensor;
144 tensor.reserve(10);
145
146 const std::vector<double> zero(num_global_cells, 0.0);
147 tensor.push_back(&zero);
148
149 static_assert(dim == 3, "");
150 std::array<int,9> kmap;
151 permeability_kind_ = fillTensor(deck, tensor, kmap);
152
153 // Assign permeability values only if such values are
154 // given in the input deck represented by 'deck'. In
155 // other words: Don't set any (arbitrary) default values.
156 // It is infinitely better to experience a reproducible
157 // crash than subtle errors resulting from a (poorly
158 // chosen) default value...
159 //
160 if (tensor.size() > 1) {
161 const int nc = global_cell.size();
162 int off = 0;
163
164 for (int c = 0; c < nc; ++c, off += dim*dim) {
165 SharedPermTensor K(dim, dim, &permeability_[off]);
166 int kix = 0;
167 const int glob = global_cell[c];
168
169 for (int i = 0; i < dim; ++i) {
170 for (int j = 0; j < dim; ++j, ++kix) {
171 K(i,j) = (*tensor[kmap[kix]])[glob];
172 }
173 K(i,i) = std::max(K(i,i), perm_threshold);
174 }
175
176 permfield_valid_[c] = std::vector<unsigned char>::value_type(1);
177 }
178 }
179 }
180
181
182
183
184
185} // namespace Opm
186
187
188#endif // OPM_ROCK_IMPL_HEADER_INCLUDED
PermTensor permeability(int cell_index) const
Read-access to permeability.
Definition: Rock_impl.hpp:84
SharedPermTensor permeabilityModifiable(int cell_index)
Read- and write-access to permeability. Use with caution.
Definition: Rock_impl.hpp:95
SharedCMatrix SharedPermTensor
Tensor type for read and write access to permeability.
Definition: Rock.hpp:43
ImmutableCMatrix PermTensor
Tensor type for read-only access to permeability.
Definition: Rock.hpp:39
Rock()
Default constructor.
Definition: Rock_impl.hpp:37
void init(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold=0.0)
Initialize from a grdecl file.
Definition: Rock_impl.hpp:44
void assignPermeability(Opm::DeckConstPtr deck, const std::vector< int > &global_cell, const double perm_threshold)
Definition: Rock_impl.hpp:132
void assignPorosity(Opm::DeckConstPtr deck, const std::vector< int > &global_cell)
Definition: Rock_impl.hpp:115
double porosity(int cell_index) const
Read-access to porosity.
Definition: Rock_impl.hpp:76
Definition: BlackoilFluid.hpp:32
@ Invalid
Definition: ReservoirPropertyCommon.hpp:50
void zero(Matrix &A)
Zero-fill a.
Definition: Matrix.hpp:603