21 #ifndef DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
22 #define DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
26 #include <opm/grid/GridHelpers.hpp>
27 #include <opm/grid/utility/OpmParserIncludes.hpp>
30 #include <opm/grid/utility/platform_dependent/disable_warnings.h>
32 #include <opm/grid/CpGrid.hpp>
33 #include <dune/common/iteratorfacades.hh>
35 #include <opm/grid/utility/platform_dependent/reenable_warnings.h>
49 : grid_(grid), cell_index_(cell_index)
54 return grid_->
faceCell(cell_index_, local_index);
86 return grid_->
faceCell(cell_index, local_index);
114 return o.index_-index_;
118 return index_==o.index_;
130 template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
136 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
141 :
IndexIterator(inner_index), grid_(grid), outer_index_(outer_index)
143 int dereference()
const
145 return std::mem_fn(AccessMethod)(*grid_, outer_index_, this->index_);
147 int elementAt(
int n)
const
149 return std::mem_fn(AccessMethod)(*grid_, outer_index_, n);
162 : grid_(grid), cell_index_(cell_index)
167 return std::mem_fn(AccessMethod)(*grid_, cell_index_, local_index);
169 const_iterator begin()
171 return const_iterator(grid_, cell_index_, 0);
175 return const_iterator(grid_, cell_index_,
176 std::mem_fn(SizeMethod)(*grid_, cell_index_));
188 template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
212 return std::mem_fn(AccessMethod)(*grid_, cell_index, local_index);
235 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
240 int index,
int cell_index)
243 int dereference()
const
245 return row_[this->index_].index();
247 int elementAt(
int n)
const
249 return row_[n].index();
251 int getCellIndex()
const
270 const int cell_index)
271 : row_(row), cell_index_(cell_index)
274 const_iterator begin()
const
276 return const_iterator(row_, 0, cell_index_);
279 const_iterator end()
const
281 return const_iterator(row_, row_.size(), cell_index_);
292 const int cell_index_;
307 auto& row=grid_->cellFaceRow(cell_index);
314 return grid_->numCellFaces();
325 namespace UgGridHelpers
333 template<const Dune::FieldVector<
double, 3>& (Dune::CpGr
id::*Method)(
int)const>
335 :
public Dune::RandomAccessIteratorFacade<CpGridCentroidIterator<Method>, Dune::FieldVector<double, 3>,
336 const Dune::FieldVector<double, 3>&,
int>
343 : grid_(&grid), cell_index_(cell_index)
346 const Dune::FieldVector<double, 3>& dereference()
const
348 return std::mem_fn(Method)(*grid_, cell_index_);
354 const Dune::FieldVector<double, 3>& elementAt(
int n)
const
356 return std::mem_fn(Method)(*grid_, n);
366 int distanceTo(
const CpGridCentroidIterator& o)
const
368 return o.cell_index_-cell_index_;
370 bool equals(
const CpGridCentroidIterator& o)
const
372 return o.grid_==grid_ && o.cell_index_==cell_index_;
384 typedef const double* ValueType;
387 typedef Dune::FieldVector<double, 3> Vector;
414 std::vector<int> createACTNUM(
const Dune::CpGrid& grid);
419 EclipseGrid createEclipseGrid(
const Dune::CpGrid& grid,
const EclipseGrid& inputGrid);
429 double cellCentroidCoordinate(
const Dune::CpGrid& grid,
int cell_index,
435 const double* cellCentroid(
const Dune::CpGrid& grid,
int cell_index);
440 double cellCenterDepth(
const Dune::CpGrid& grid,
int cell_index);
456 Vector faceAreaNormalEcl(
const Dune::CpGrid& grid,
int face_index);
461 double cellVolume(
const Dune::CpGrid& grid,
int cell_index);
465 :
public Dune::RandomAccessIteratorFacade<CellVolumeIterator, double, double, int>
472 : grid_(&grid), cell_index_(cell_index)
475 double dereference()
const
477 return grid_->cellVolume(cell_index_);
483 double elementAt(
int n)
const
485 return grid_->cellVolume(n);
495 int distanceTo(
const CellVolumeIterator& o)
const
497 return o.cell_index_-cell_index_;
499 bool equals(
const CellVolumeIterator& o)
const
501 return o.grid_==grid_ && o.cell_index_==cell_index_;
525 typedef const Dune::CpGrid::Vector ValueType;
564 const double* vertexCoordinates(
const Dune::CpGrid& grid,
int index);
566 const double* faceNormal(
const Dune::CpGrid& grid,
int face_index);
568 double faceArea(
const Dune::CpGrid& grid,
int face_index);
[ provides Dune::Grid ]
Definition: CpGrid.hpp:207
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition: CpGrid.hpp:995
Definition: GridHelpers.hpp:297
std::size_t noEntries() const
Get the number of non-zero entries.
Definition: GridHelpers.hpp:312
Definition: GridHelpers.hpp:237
Definition: GridHelpers.hpp:232
A class representing the face to cells mapping similar to the way done in UnstructuredGrid.
Definition: GridHelpers.hpp:64
FaceCellsContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:70
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition: GridHelpers.hpp:84
FaceCellsProxy operator[](int cell_index) const
Get the mapping for a cell.
Definition: GridHelpers.hpp:75
A proxy class representing a row of FaceCellsContainer.
Definition: GridHelpers.hpp:43
FaceCellsProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition: GridHelpers.hpp:48
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition: GridHelpers.hpp:52
A class representing the face to vertices mapping similar to the way done in UnstructuredGrid.
Definition: GridHelpers.hpp:222
FaceVerticesContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:226
Definition: GridHelpers.hpp:94
A class representing the sparse mapping of entity relations (e.g.
Definition: GridHelpers.hpp:191
LocalIndexContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition: GridHelpers.hpp:196
row_type operator[](int cell_index) const
Get the mapping for a cell.
Definition: GridHelpers.hpp:201
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition: GridHelpers.hpp:210
Definition: GridHelpers.hpp:138
A proxy class representing a row of LocalIndexContainerProxy.
Definition: GridHelpers.hpp:133
LocalIndexProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition: GridHelpers.hpp:161
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition: GridHelpers.hpp:165
A class used as a row type for OrientedEntityTable.
Definition: OrientedEntityTable.hpp:55
An iterator over the cell volumes.
Definition: GridHelpers.hpp:466
CellVolumeIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition: GridHelpers.hpp:471
An iterator over the cell volumes.
Definition: GridHelpers.hpp:337
CpGridCentroidIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition: GridHelpers.hpp:342
Copyright 2019 Equinor AS.
Definition: CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:29
face_tag
Connection taxonomy.
Definition: preprocess.h:66
Definition: GridHelpers.hpp:329
Traits of the cell centroids of a grid.
Definition: GridHelpers.hpp:126
The mapping of the grid type to type of the iterator over the cell volumes.
Definition: GridHelpers.hpp:191
Maps the grid type to the associated type of the face to vertices mapping.
Definition: GridHelpers.hpp:292
Traits of the face to attached cell mappping of a grid.
Definition: GridHelpers.hpp:334
Traits of the face centroids of a grid.
Definition: GridHelpers.hpp:241