21 #ifndef DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
22 #define DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
25 #include <opm/core/grid/GridHelpers.hpp>
26 #include <opm/grid/utility/OpmParserIncludes.hpp>
29 #include <opm/grid/utility/platform_dependent/disable_warnings.h>
31 #include <dune/grid/CpGrid.hpp>
32 #include <dune/common/iteratorfacades.hh>
34 #include <opm/grid/utility/platform_dependent/reenable_warnings.h>
48 : grid_(grid), cell_index_(cell_index)
53 return grid_->
faceCell(cell_index_, local_index);
85 return grid_->
faceCell(cell_index, local_index);
113 return o.index_-index_;
117 return index_==o.index_;
129 template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
135 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
140 :
IndexIterator(inner_index), grid_(grid), outer_index_(outer_index)
142 int dereference()
const
144 return std::mem_fn(AccessMethod)(*grid_, outer_index_, this->index_);
146 int elementAt(
int n)
const
148 return std::mem_fn(AccessMethod)(*grid_, outer_index_, n);
161 : grid_(grid), cell_index_(cell_index)
166 return std::mem_fn(AccessMethod)(*grid_, cell_index_, local_index);
168 const_iterator begin()
170 return const_iterator(grid_, cell_index_, 0);
174 return const_iterator(grid_, cell_index_,
175 std::mem_fn(SizeMethod)(*grid_, cell_index_));
187 template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
202 return row_type(grid_, cell_index);
211 return std::mem_fn(AccessMethod)(*grid_, cell_index, local_index);
234 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
239 int index,
int cell_index)
242 int dereference()
const
244 return row_[this->index_].index();
246 int elementAt(
int n)
const
248 return row_[n].index();
250 int getCellIndex()
const
269 const int cell_index)
270 : row_(row), cell_index_(cell_index)
273 const_iterator begin()
const
275 return const_iterator(row_, 0, cell_index_);
278 const_iterator end()
const
280 return const_iterator(row_, row_.size(), cell_index_);
291 const int cell_index_;
306 auto& row=grid_->cellFaceRow(cell_index);
313 return grid_->numCellFaces();
324 namespace UgGridHelpers
332 template<const Dune::FieldVector<
double, 3>& (Dune::CpGr
id::*Method)(
int)const>
334 :
public Dune::RandomAccessIteratorFacade<CpGridCentroidIterator<Method>, Dune::FieldVector<double, 3>,
335 const Dune::FieldVector<double, 3>&,
int>
342 : grid_(&grid), cell_index_(cell_index)
345 const Dune::FieldVector<double, 3>& dereference()
const
347 return std::mem_fn(Method)(*grid_, cell_index_);
353 const Dune::FieldVector<double, 3>& elementAt(
int n)
const
355 return std::mem_fn(Method)(*grid_, n);
365 int distanceTo(
const CpGridCentroidIterator& o)
const
367 return o.cell_index_-cell_index_;
369 bool equals(
const CpGridCentroidIterator& o)
const
371 return o.grid_==grid_ && o.cell_index_==cell_index_;
383 typedef const double* ValueType;
386 typedef Dune::FieldVector<double, 3> Vector;
410 EclipseGrid createEclipseGrid(
const Dune::CpGrid& grid,
const EclipseGrid& inputGrid);
422 double cellCentroidCoordinate(
const Dune::CpGrid& grid,
int cell_index,
428 const double* cellCentroid(
const Dune::CpGrid& grid,
int cell_index);
433 double cellCenterDepth(
const Dune::CpGrid& grid,
int cell_index);
449 Vector faceAreaNormalEcl(
const Dune::CpGrid& grid,
int face_index);
454 double cellVolume(
const Dune::CpGrid& grid,
int cell_index);
458 :
public Dune::RandomAccessIteratorFacade<CellVolumeIterator, double, double, int>
465 : grid_(&grid), cell_index_(cell_index)
468 double dereference()
const
470 return grid_->cellVolume(cell_index_);
476 double elementAt(
int n)
const
478 return grid_->cellVolume(n);
488 int distanceTo(
const CellVolumeIterator& o)
const
490 return o.cell_index_-cell_index_;
492 bool equals(
const CellVolumeIterator& o)
const
494 return o.grid_==grid_ && o.cell_index_==cell_index_;
518 typedef const Dune::CpGrid::Vector ValueType;
557 const double* vertexCoordinates(
const Dune::CpGrid& grid,
int index);
559 const double* faceNormal(
const Dune::CpGrid& grid,
int face_index);
561 double faceArea(
const Dune::CpGrid& grid,
int face_index);
FaceCellsProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition: GridHelpers.hpp:47
CpGridCentroidIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition: GridHelpers.hpp:341
row_type operator[](int cell_index) const
Get the mapping for a cell.
Definition: GridHelpers.hpp:200
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition: CpGrid.hpp:769
FaceVerticesContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:225
A proxy class representing a row of FaceCellsContainer.
Definition: GridHelpers.hpp:41
FaceCellsContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:69
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition: GridHelpers.hpp:209
A proxy class representing a row of LocalIndexContainerProxy.
Definition: GridHelpers.hpp:131
Definition: GridHelpers.hpp:233
[ provides Dune::Grid ]
Definition: CpGrid.hpp:213
CellVolumeIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition: GridHelpers.hpp:464
A class representing the sparse mapping of entity relations (e.g.
Definition: GridHelpers.hpp:189
An iterator over the cell volumes.
Definition: GridHelpers.hpp:457
Definition: GridHelpers.hpp:327
Definition: GridHelpers.hpp:295
FaceCellsProxy operator[](int cell_index) const
Get the mapping for a cell.
Definition: GridHelpers.hpp:74
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition: GridHelpers.hpp:83
Definition: GridHelpers.hpp:230
Traits of the cell centroids of a grid.
Definition: GridHelpers.hpp:131
A class representing the face to cells mapping similar to the way done in UnstructuredGrid.
Definition: GridHelpers.hpp:62
LocalIndexProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition: GridHelpers.hpp:160
Traits of the face to attached cell mappping of a grid.
Definition: GridHelpers.hpp:345
std::size_t noEntries() const
Get the number of non-zero entries.
Definition: GridHelpers.hpp:311
A class representing the face to vertices mapping similar to the way done in UnstructuredGrid.
Definition: GridHelpers.hpp:219
An iterator over the cell volumes.
Definition: GridHelpers.hpp:333
Maps the grid type to the associated type of the face to vertices mapping.
Definition: GridHelpers.hpp:303
face_tag
Connection taxonomy.
Definition: preprocess.h:67
Definition: GridHelpers.hpp:92
Traits of the face centroids of a grid.
Definition: GridHelpers.hpp:246
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition: GridHelpers.hpp:51
LocalIndexContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:195
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
Definition: GridHelpers.hpp:134
The mapping of the grid type to type of the iterator over the cell volumes.
Definition: GridHelpers.hpp:196
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition: GridHelpers.hpp:164