36 #ifndef OPM_ENTITY_HEADER
37 #define OPM_ENTITY_HEADER
39 #include <dune/common/version.hh>
40 #include <dune/geometry/type.hh>
41 #include <dune/grid/common/gridenums.hh>
43 #include "PartitionTypeIndicator.hpp"
44 #include "EntityRep.hpp"
51 template<
int,
int>
class Geometry;
52 template<
int,PartitionIteratorType>
class Iterator;
53 class IntersectionIterator;
54 class HierarchicIterator;
56 class LevelGlobalIdSet;
71 enum { codimension = codim };
72 enum { dimension = 3 };
73 enum { mydimension = dimension - codimension };
74 enum { dimensionworld = 3 };
115 :
EntityRep<codim>(entityrep), pgrid_(&grid)
121 :
EntityRep<codim>(index_arg, orientation_arg), pgrid_(&grid)
171 return Dune::GeometryTypes::cube(3 - codim);
243 const Entity& impl()
const
266 #include "Iterators.hpp"
267 #include "Entity.hpp"
268 #include "Intersection.hpp"
276 static_assert(codim == 0,
"");
283 static_assert(codim == 0,
"");
290 static_assert(codim == 0,
"");
297 static_assert(codim == 0,
"");
318 return pgrid_->partition_type_indicator_->getPartitionType(*
this);
323 #include <opm/grid/cpgrid/CpGridData.hpp>
324 #include <opm/grid/cpgrid/Intersection.hpp>
331 inline unsigned int numFaces(
const OrientedEntityTable<0, 1>& cell_to_face,
334 return cell_to_face[e].size();
338 unsigned int numFaces(
const OrientedEntityTable<0, 1>&,
const Entity<codim>&)
350 }
else if ( codim == 0 ){
352 return Detail::numFaces(pgrid_->cell_to_face_, *
this);
353 }
else if ( cc == 3 ) {
363 return pgrid_->geomVector<codim>()[*
this];
370 static_assert(codim == 0,
"");
375 }
else if (cc == 3) {
376 assert(i >= 0 && i < 8);
378 typename Codim<cc>::Entity se(*pgrid_, corner_index,
true);
381 OPM_THROW(std::runtime_error,
"No subentity exists of codimension " << cc);
391 Iter end = ileafend();
392 for (Iter it = ileafbegin(); it != end; ++it) {
393 if (it->boundary())
return true;
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:123
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:98
bool operator==(const EntityRep &other) const
Equality operator.
Definition: EntityRep.hpp:178
int index() const
The (positive) index of an entity.
Definition: EntityRep.hpp:125
Definition: Entity.hpp:64
bool mightVanish() const
Returns true, if entity might disappear during the next call to adapt(). Dummy.
Definition: Entity.hpp:210
Codim< cc >::Entity subEntity(int i) const
Obtain subentity.
Entity father() const
Dummy, returning this.
Definition: Entity.hpp:223
unsigned int subEntities(const unsigned int cc) const
The count of subentities of codimension cc.
Definition: Entity.hpp:345
LeafIntersectionIterator ileafend() const
End iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:295
LocalGeometry geometryInFather() const
Dummy, returning default geometry.
Definition: Entity.hpp:230
EntitySeed seed() const
Return an entity seed.
Definition: Entity.hpp:139
HierarchicIterator hend(int) const
Dummy beyond last child iterator.
Definition: Entity.hpp:310
const Geometry & geometry() const
Returns the geometry of the entity (does not depend on its orientation).
Definition: Entity.hpp:361
Entity(const CpGridData &grid, int index_arg, bool orientation_arg)
Constructor taking a grid, entity index and orientation.
Definition: Entity.hpp:120
bool hasFather() const
No hierarchy in this grid.
Definition: Entity.hpp:216
bool isValid() const
isValid method for EntitySeed
Definition: Entity.hpp:399
bool operator==(const Entity &other) const
Equality.
Definition: Entity.hpp:126
bool isRegular() const
Refinement is not defined for CpGrid.
Definition: Entity.hpp:160
int count() const
The count of subentities of codimension cc.
Definition: Entity.hpp:179
Entity()
Constructor taking a grid and an integer entity representation.
Definition: Entity.hpp:108
HierarchicIterator hbegin(int) const
Dummy first child iterator.
Definition: Entity.hpp:303
PartitionType partitionType() const
For now, the grid is serial and the only partitionType() is InteriorEntity.
Definition: Entity.hpp:316
LeafIntersectionIterator ileafbegin() const
Start iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:288
Entity(const CpGridData &grid, EntityRep< codim > entityrep)
Constructor taking a grid and an entity representation.
Definition: Entity.hpp:114
LevelIntersectionIterator ilevelend() const
End iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:281
bool isNew() const
Returns true, if the entity has been created during the last call to adapt(). Dummy.
Definition: Entity.hpp:204
LevelIntersectionIterator ilevelbegin() const
Start iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:274
GeometryType type() const
Using the cube type for all entities now (cells and vertices).
Definition: Entity.hpp:169
int level() const
We do not support refinement, so level() is always 0.
Definition: Entity.hpp:148
bool isLeaf() const
The entity is always on the leaf grid, since we have no refinement.
Definition: Entity.hpp:154
bool operator!=(const Entity &other) const
Inequality.
Definition: Entity.hpp:132
bool hasBoundaryIntersections() const
Returns true if any of my intersections are on the boundary.
Definition: Entity.hpp:386
This class encapsulates geometry for both vertices, intersections and cells.
Definition: Geometry.hpp:71
The global id set for Dune.
Definition: Indexsets.hpp:325
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:104
Definition: Intersection.hpp:288
Definition: Indexsets.hpp:257
Copyright 2019 Equinor AS.
Definition: CartesianIndexMapper.hpp:10
Definition: Entity.hpp:84