21#ifndef DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
22#define DUNE_CORNERPOINT_GRIDHELPERS_HEADER_INCLUDED
26#include <opm/grid/GridHelpers.hpp>
28#include <opm/grid/utility/platform_dependent/disable_warnings.h>
30#include <opm/grid/CpGrid.hpp>
31#include <dune/common/iteratorfacades.hh>
33#include <opm/grid/utility/platform_dependent/reenable_warnings.h>
47 : grid_(grid), cell_index_(cell_index)
52 return grid_->
faceCell(cell_index_, local_index);
84 return grid_->
faceCell(cell_index, local_index);
112 return o.index_-index_;
116 return index_==o.index_;
128template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
134 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
139 :
IndexIterator(inner_index), grid_(grid), outer_index_(outer_index)
141 int dereference()
const
143 return std::mem_fn(AccessMethod)(*grid_, outer_index_, this->index_);
145 int elementAt(
int n)
const
147 return std::mem_fn(AccessMethod)(*grid_, outer_index_, n);
160 : grid_(grid), cell_index_(cell_index)
165 return std::mem_fn(AccessMethod)(*grid_, cell_index_, local_index);
167 const_iterator begin()
169 return const_iterator(grid_, cell_index_, 0);
173 return const_iterator(grid_, cell_index_,
174 std::mem_fn(SizeMethod)(*grid_, cell_index_));
186template<int (
Dune::CpGrid::*AccessMethod)(int,int)
const,
210 return std::mem_fn(AccessMethod)(*grid_, cell_index, local_index);
233 :
public Dune::RandomAccessIteratorFacade<iterator,int, int, int>,
238 int index,
int cell_index)
241 int dereference()
const
243 return row_[this->index_].index();
245 int elementAt(
int n)
const
247 return row_[n].index();
249 int getCellIndex()
const
268 const int cell_index)
269 : row_(row), cell_index_(cell_index)
272 const_iterator begin()
const
274 return const_iterator(row_, 0, cell_index_);
277 const_iterator end()
const
279 return const_iterator(row_, row_.size(), cell_index_);
290 const int cell_index_;
305 auto& row=grid_->cellFaceRow(cell_index);
312 return grid_->numCellFaces();
323namespace UgGridHelpers
331template<const Dune::FieldVector<
double, 3>& (Dune::CpGr
id::*Method)(
int)const>
333 :
public Dune::RandomAccessIteratorFacade<CpGridCentroidIterator<Method>, Dune::FieldVector<double, 3>,
334 const Dune::FieldVector<double, 3>&,
int>
341 : grid_(&grid), cell_index_(cell_index)
344 const Dune::FieldVector<double, 3>& dereference()
const
346 return std::mem_fn(Method)(*grid_, cell_index_);
352 const Dune::FieldVector<double, 3>& elementAt(
int n)
const
354 return std::mem_fn(Method)(*grid_, n);
364 int distanceTo(
const CpGridCentroidIterator& o)
const
366 return o.cell_index_-cell_index_;
368 bool equals(
const CpGridCentroidIterator& o)
const
370 return o.grid_==grid_ && o.cell_index_==cell_index_;
382 typedef const double* ValueType;
385typedef Dune::FieldVector<double, 3> Vector;
417EclipseGrid createEclipseGrid(
const Dune::CpGrid& grid,
const EclipseGrid& inputGrid);
427double cellCentroidCoordinate(
const Dune::CpGrid& grid,
int cell_index,
433const double* cellCentroid(
const Dune::CpGrid& grid,
int cell_index);
438double cellCenterDepth(
const Dune::CpGrid& grid,
int cell_index);
454Vector faceAreaNormalEcl(
const Dune::CpGrid& grid,
int face_index);
459double cellVolume(
const Dune::CpGrid& grid,
int cell_index);
463 :
public Dune::RandomAccessIteratorFacade<CellVolumeIterator, double, double, int>
470 : grid_(&grid), cell_index_(cell_index)
473 double dereference()
const
475 return grid_->cellVolume(cell_index_);
481 double elementAt(
int n)
const
483 return grid_->cellVolume(n);
493 int distanceTo(
const CellVolumeIterator& o)
const
495 return o.cell_index_-cell_index_;
497 bool equals(
const CellVolumeIterator& o)
const
499 return o.grid_==grid_ && o.cell_index_==cell_index_;
523 typedef const Dune::CpGrid::Vector ValueType;
562const double* vertexCoordinates(
const Dune::CpGrid& grid,
int index);
564const double* faceNormal(
const Dune::CpGrid& grid,
int face_index);
566double faceArea(
const Dune::CpGrid& grid,
int face_index);
[ provides Dune::Grid ]
Definition CpGrid.hpp:201
int faceCell(int face, int local_index) const
Get the index identifying a cell attached to a face.
Definition CpGrid.cpp:1113
Definition GridHelpers.hpp:295
std::size_t noEntries() const
Get the number of non-zero entries.
Definition GridHelpers.hpp:310
Definition GridHelpers.hpp:235
Definition GridHelpers.hpp:230
A class representing the face to cells mapping similar to the way done in UnstructuredGrid.
Definition GridHelpers.hpp:62
FaceCellsContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition GridHelpers.hpp:68
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition GridHelpers.hpp:82
FaceCellsProxy operator[](int cell_index) const
Get the mapping for a cell.
Definition GridHelpers.hpp:73
A proxy class representing a row of FaceCellsContainer.
Definition GridHelpers.hpp:41
FaceCellsProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition GridHelpers.hpp:46
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition GridHelpers.hpp:50
A class representing the face to vertices mapping similar to the way done in UnstructuredGrid.
Definition GridHelpers.hpp:220
FaceVerticesContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition GridHelpers.hpp:224
Definition GridHelpers.hpp:92
A class representing the sparse mapping of entity relations (e.g.
Definition GridHelpers.hpp:189
LocalIndexContainerProxy(const Dune::CpGrid *grid)
Constructor.
Definition GridHelpers.hpp:194
row_type operator[](int cell_index) const
Get the mapping for a cell.
Definition GridHelpers.hpp:199
int operator()(int cell_index, int local_index) const
Get a face associated with a cell.
Definition GridHelpers.hpp:208
Definition GridHelpers.hpp:136
A proxy class representing a row of LocalIndexContainerProxy.
Definition GridHelpers.hpp:131
LocalIndexProxy(const Dune::CpGrid *grid, int cell_index)
Constructor.
Definition GridHelpers.hpp:159
int operator[](int local_index)
Get the index of the cell associated with a local_index.
Definition GridHelpers.hpp:163
A class used as a row type for OrientedEntityTable.
Definition OrientedEntityTable.hpp:55
An iterator over the cell volumes.
Definition GridHelpers.hpp:464
CellVolumeIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition GridHelpers.hpp:469
An iterator over the cell volumes.
Definition GridHelpers.hpp:335
CpGridCentroidIterator(const Dune::CpGrid &grid, int cell_index)
Creates an iterator.
Definition GridHelpers.hpp:340
Copyright 2019 Equinor AS.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68
face_tag
Connection taxonomy.
Definition preprocess.h:66
Maps the grid type to the associated type of the cell to faces mapping.
Definition GridHelpers.hpp:280
Traits of the cell centroids of a grid.
Definition GridHelpers.hpp:129
The mapping of the grid type to type of the iterator over the cell volumes.
Definition GridHelpers.hpp:194
Maps the grid type to the associated type of the face to vertices mapping.
Definition GridHelpers.hpp:295
Traits of the face to attached cell mappping of a grid.
Definition GridHelpers.hpp:337
Traits of the face centroids of a grid.
Definition GridHelpers.hpp:244