GridHelpers.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2014, 2015 Dr. Markus Blatt - HPC-Simulation-Software & Services
3  Copyright 2014 Statoil AS
4  Copyright 2015
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 3 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
21 #ifndef OPM_CORE_GRIDHELPERS_HEADER_INCLUDED
22 #define OPM_CORE_GRIDHELPERS_HEADER_INCLUDED
23 #include <opm/core/grid.h>
24 #include <boost/range/iterator_range.hpp>
25 
26 namespace Opm
27 {
28 namespace UgGridHelpers
29 {
30 
37 {
38 public:
39  class IntRange : public boost::iterator_range<const int*>
40  {
41  public:
42  typedef boost::iterator_range<const int*> BaseRowType;
43  typedef BaseRowType::size_type size_type;
44  typedef int value_type;
45 
46  IntRange(const int* start_arg, const int* end_arg)
47  : BaseRowType(start_arg, end_arg)
48  {}
49  };
51  typedef boost::iterator_range<const int*> row_type;
52 
58  SparseTableView(int* data, int *offset, std::size_t size_arg)
59  : data_(data), offset_(offset), size_(size_arg)
60  {}
61 
65  row_type operator[](std::size_t row) const
66  {
67  assert(row<=size());
68  return row_type(data_ + offset_[row], data_ + offset_[row+1]);
69  }
70 
73  std::size_t size() const
74  {
75  return size_;
76  }
77 
79  std::size_t noEntries() const
80  {
81  return offset_[size_];
82  }
83 
84 private:
86  const int* data_;
90  const int* offset_;
92  std::size_t size_;
93 };
94 
96 int numCells(const UnstructuredGrid& grid);
97 
99 int numFaces(const UnstructuredGrid& grid);
100 
102 int dimensions(const UnstructuredGrid& grid);
103 
106 
108 const int* cartDims(const UnstructuredGrid& grid);
109 
114 const int* globalCell(const UnstructuredGrid& grid);
115 
121 template<class G>
123 {
124 };
125 
126 template<>
128 {
129  typedef const double* IteratorType;
130  typedef const double* ValueType;
131 };
132 
139 
140 
144 double cellCenterDepth(const UnstructuredGrid& grid, int cell_index);
145 
146 
151 double cellCentroidCoordinate(const UnstructuredGrid& grid, int cell_index,
152  int coordinate);
153 
154 
158 const double* cellCentroid(const UnstructuredGrid& grid, int cell_index);
159 
160 
164 double cellVolume(const UnstructuredGrid& grid, int cell_index);
165 
171 template<class T>
173 {
174 };
175 
176 template<>
178 {
179  typedef const double* IteratorType;
180 };
181 
183 const double* beginCellVolumes(const UnstructuredGrid& grid);
184 
186 const double* endCellVolumes(const UnstructuredGrid& grid);
187 
191 const double* faceCentroid(const UnstructuredGrid& grid, int face_index);
192 
193 
199 template<class G>
201 {
202 };
203 
204 template<>
206 {
207  typedef const double* IteratorType;
208  typedef const double* ValueType;
209 };
210 
214 
220 faceCentroid(const UnstructuredGrid& grid, int face_index);
221 
225 const double* faceNormal(const UnstructuredGrid& grid, int face_index);
226 
230 double faceArea(const UnstructuredGrid& grid, int face_index);
231 
236 int faceTag(const UnstructuredGrid& grid, boost::iterator_range<const int*>::const_iterator cell_face);
237 
242 template<class T>
244 {
245 };
246 
247 template<>
249 {
251 };
252 
257 template<class T>
259 {
260 };
261 
262 template<>
264 {
266 };
267 
271 
275 
279 const double* vertexCoordinates(const UnstructuredGrid& grid, int index);
280 
282 {
283 public:
285  : face_cells_(grid.face_cells)
286  {}
287  int operator()(int face_index, int local_index) const
288  {
289  return face_cells_[2*face_index+local_index];
290  }
291 private:
292  const int* face_cells_;
293 };
294 
299 template<class T>
301 {};
302 
303 template<>
305 {
307 };
308 
311 
317 template<class T>
318 T* increment(T* cc, int i, int dim)
319 {
320  return cc+(i*dim);
321 }
326 template<class T>
327 T increment(const T& t, int i, int)
328 {
329  return t+i;
330 }
331 
336 template<class T>
337 double getCoordinate(T* cc, int i)
338 {
339  return cc[i];
340 }
341 
347 template<class T>
348 double getCoordinate(T t, int i)
349 {
350  return (*t)[i];
351 }
352 
353 } // end namespace UGGridHelpers
354 } // end namespace OPM
355 #endif
The mapping of the grid type to type of the iterator over the cell volumes.
Definition: GridHelpers.hpp:172
BaseRowType::size_type size_type
Definition: GridHelpers.hpp:43
const double * ValueType
Definition: GridHelpers.hpp:130
IntRange(const int *start_arg, const int *end_arg)
Definition: GridHelpers.hpp:46
std::size_t noEntries() const
Get the number of non-zero entries.
Definition: GridHelpers.hpp:79
boost::iterator_range< const int * > row_type
The type of the roww.
Definition: GridHelpers.hpp:51
const double * beginCellVolumes(const UnstructuredGrid &grid)
Get an iterator over the cell volumes of a grid positioned at the first cell.
Definition: grid.h:98
double getCoordinate(T *cc, int i)
Get the i-th corrdinate of a centroid.
Definition: GridHelpers.hpp:337
FaceCellTraits< UnstructuredGrid >::Type faceCells(const UnstructuredGrid &grid)
Get the face to cell mapping of a grid.
Definition: AnisotropicEikonal.hpp:43
Definition: GridHelpers.hpp:281
FaceCellsProxy Type
Definition: GridHelpers.hpp:306
Face2VerticesTraits< UnstructuredGrid >::Type face2Vertices(const UnstructuredGrid &grid)
Get the face to vertices mapping of a grid.
FaceCellsProxy(const UnstructuredGrid &grid)
Definition: GridHelpers.hpp:284
int value_type
Definition: GridHelpers.hpp:44
int faceTag(const UnstructuredGrid &grid, boost::iterator_range< const int * >::const_iterator cell_face)
Get Eclipse Cartesian tag of a face.
Maps the grid type to the associated type of the face to vertices mapping.
Definition: GridHelpers.hpp:258
Traits of the face centroids of a grid.
Definition: GridHelpers.hpp:200
const double * cellCentroid(const UnstructuredGrid &grid, int cell_index)
Get the centroid of a cell.
const double * faceNormal(const UnstructuredGrid &grid, int face_index)
Get the normal of a face.
int numCells(const UnstructuredGrid &grid)
Get the number of cells of a grid.
int numFaces(const UnstructuredGrid &grid)
Get the number of faces of a grid.
row_type operator[](std::size_t row) const
Get a row of the the table.
Definition: GridHelpers.hpp:65
Traits of the cell centroids of a grid.
Definition: GridHelpers.hpp:122
const int * globalCell(const UnstructuredGrid &grid)
Get the local to global index mapping.
CellCentroidTraits< UnstructuredGrid >::IteratorType beginCellCentroids(const UnstructuredGrid &grid)
Get an iterator over the cell centroids positioned at the first cell.
const int * cartDims(const UnstructuredGrid &grid)
Get the cartesion dimension of the underlying structured grid.
double cellCentroidCoordinate(const UnstructuredGrid &grid, int cell_index, int coordinate)
Get a coordinate of a specific cell centroid.
int dimensions(const UnstructuredGrid &grid)
Get the dimensions of a grid.
const double * IteratorType
Definition: GridHelpers.hpp:179
double cellVolume(const UnstructuredGrid &grid, int cell_index)
Get the volume of a cell.
const double * ValueType
Definition: GridHelpers.hpp:208
boost::iterator_range< const int * > BaseRowType
Definition: GridHelpers.hpp:42
int operator()(int face_index, int local_index) const
Definition: GridHelpers.hpp:287
Maps the grid type to the associated type of the cell to faces mapping.
Definition: GridHelpers.hpp:243
T * increment(T *cc, int i, int dim)
Increment an iterator over an array that reresents a dense row-major matrix with dims columns...
Definition: GridHelpers.hpp:318
const double * IteratorType
Definition: GridHelpers.hpp:129
std::size_t size() const
Get the size of the table.
Definition: GridHelpers.hpp:73
FaceCentroidTraits< UnstructuredGrid >::IteratorType beginFaceCentroids(const UnstructuredGrid &grid)
Get an iterator over the face centroids positioned at the first cell.
const UnstructuredGrid & grid
Definition: ColumnExtract.hpp:31
const double * IteratorType
Definition: GridHelpers.hpp:207
double cellCenterDepth(const UnstructuredGrid &grid, int cell_index)
Get vertical position of cell center ("zcorn" average.)
Definition: GridHelpers.hpp:39
SparseTableView(int *data, int *offset, std::size_t size_arg)
Creates a sparse table view.
Definition: GridHelpers.hpp:58
SparseTableView Type
Definition: GridHelpers.hpp:250
Cell2FacesTraits< UnstructuredGrid >::Type cell2Faces(const UnstructuredGrid &grid)
Get the cell to faces mapping of a grid.
Allows viewing a sparse table consisting out of C-array.
Definition: GridHelpers.hpp:36
const double * faceCentroid(const UnstructuredGrid &grid, int face_index)
Get the cell centroid of a face.
Traits of the face to attached cell mappping of a grid.
Definition: GridHelpers.hpp:300
const double * endCellVolumes(const UnstructuredGrid &grid)
Get an iterator over the cell volumes of a grid positioned after the last cell.
const double * vertexCoordinates(const UnstructuredGrid &grid, int index)
Get the coordinates of a vertex of the grid.
double faceArea(const UnstructuredGrid &grid, int face_index)
Get the area of a face.
int numCellFaces(const UnstructuredGrid &grid)
Get the number of faces, where each face counts as many times as there are adjacent faces...
SparseTableView Type
Definition: GridHelpers.hpp:265