37#ifndef OPM_ENTITY_HEADER
38#define OPM_ENTITY_HEADER
40#include <dune/common/version.hh>
41#include <dune/geometry/type.hh>
42#include <dune/geometry/referenceelements.hh>
43#include <dune/grid/common/gridenums.hh>
45#include "PartitionTypeIndicator.hpp"
46#include <opm/grid/cpgrid/DefaultGeometryPolicy.hpp>
50 const std::vector<std::array<int,3>>&,
51 const std::vector<std::array<int,3>>&,
52 const std::vector<std::array<int,3>>&,
53 const std::vector<std::string>&);
61template<
int,PartitionIteratorType>
class Iterator;
73 friend class LevelGlobalIdSet;
74 friend class GlobalIdSet;
75 friend class HierarchicIterator;
76 friend class CpGridData;
78 const std::vector<std::array<int,3>>&,
79 const std::vector<std::array<int,3>>&,
80 const std::vector<std::array<int,3>>&,
81 const std::vector<std::string>&);
87 static constexpr int dimension = 3;
88 static constexpr int mydimension = dimension -
codimension;
89 static constexpr int dimensionworld = 3;
103 typedef Geometry LocalGeometry;
109 typedef double ctype;
128 :
EntityRep<codim>(entityrep), pgrid_(&grid)
133 Entity(
const CpGridData& grid,
int index_arg,
bool orientation_arg)
134 :
EntityRep<codim>(index_arg, orientation_arg), pgrid_(&grid)
139 Entity(
int index_arg,
bool orientation_arg)
140 :
EntityRep<codim>(index_arg, orientation_arg), pgrid_()
161 return EntitySeed(
impl() );
194 return Dune::GeometryTypes::cube(3 - codim);
222 HierarchicIterator
hend(
int)
const;
295 int getIdxInParentCell()
const;
298 const CpGridData* pgrid_;
305 return Dune::referenceElement<double,3-codim>(Dune::GeometryTypes::cube(3-codim));
314#include "Iterators.hpp"
315#include "Intersection.hpp"
323 static_assert(codim == 0,
"");
324 return LevelIntersectionIterator(*pgrid_, *
this,
false);
330 static_assert(codim == 0,
"");
331 return LevelIntersectionIterator(*pgrid_, *
this,
true);
337 static_assert(codim == 0,
"");
338 return LeafIntersectionIterator(*pgrid_, *
this,
false);
344 static_assert(codim == 0,
"");
345 return LeafIntersectionIterator(*pgrid_, *
this,
true);
353 return HierarchicIterator(*
this, maxLevel);
361 return HierarchicIterator(maxLevel);
367 return pgrid_->partition_type_indicator_->getPartitionType(*
this);
374#include <opm/grid/cpgrid/CpGridData.hpp>
385 else if ( codim == 0 ){
396 return pgrid_->geomVector<codim>()[*
this];
403 static_assert(codim == 0,
"");
406 typename Codim<cc>::Entity se(*pgrid_, this->index(), this->orientation());
408 }
else if (cc == 3) {
409 assert(i >= 0 && i < 8);
410 int corner_index = pgrid_->cell_to_point_[this->index()][i];
411 typename Codim<cc>::Entity se(*pgrid_, corner_index,
true);
415 OPM_THROW(std::runtime_error,
416 "No subentity exists of codimension " + std::to_string(cc));
425 typedef LeafIntersectionIterator Iter;
427 for (Iter it =
ileafbegin(); it != end; ++it) {
428 if (it->boundary())
return true;
454 return pgrid_->leaf_to_level_cells_.empty()? pgrid_->level_ : pgrid_->leaf_to_level_cells_[this->
index()][0];
468 if (pgrid_ -> parent_to_children_cells_.empty()){
472 return (std::get<0>((pgrid_ -> parent_to_children_cells_)[this->
index()]) == -1);
479 int data_count = pgrid_->level_data_ptr_->size();
480 if (data_count == 1) {
484 assert(data_count >= 2);
489 return isLeaf() && (pgrid_->refinement_max_level_ == data_count - 2);
502 if ((pgrid_ -> child_to_parent_cells_.empty()) || (pgrid_ -> child_to_parent_cells_[this->index()][0] == -1)){
514 const int& coarser_level = pgrid_ -> child_to_parent_cells_[this->
index()][0];
515 const int& parent_cell_index = pgrid_ -> child_to_parent_cells_[this->
index()][1];
516 return Entity<0>( *((*(pgrid_ -> level_data_ptr_))[coarser_level].get()), parent_cell_index,
true);
519 OPM_THROW(std::logic_error,
"Entity has no father.");
524int Dune::cpgrid::Entity<codim>::getIdxInParentCell()
const
526 return pgrid_ -> cell_to_idxInParentCell_[this->index()];
534 OPM_THROW(std::logic_error,
"Entity has no father.");
538 static constexpr std::array<int,8> in_father_reference_elem_corner_indices = {0,1,2,3,4,5,6,7};
542 auto idx_in_parent_cell = pgrid_ -> cell_to_idxInParentCell_[this->
index()];
543 if (idx_in_parent_cell !=-1) {
544 const auto& cells_per_dim = (*(pgrid_ -> level_data_ptr_))[this->
level()] -> cells_per_dim_;
545 const auto& auxArr = pgrid_ -> getReferenceRefinedCorners(idx_in_parent_cell, cells_per_dim);
546 FieldVector<double, 3> corners_in_father_reference_elem_temp[8] =
547 { auxArr[0], auxArr[1], auxArr[2], auxArr[3], auxArr[4], auxArr[5], auxArr[6], auxArr[7]};
548 auto in_father_reference_elem_corners = std::make_shared<EntityVariable<cpgrid::Geometry<0, 3>, 3>>();
551 mutable_in_father_reference_elem_corners.assign(corners_in_father_reference_elem_temp,
552 corners_in_father_reference_elem_temp + 8);
554 Dune::FieldVector<double, 3> center_in_father_reference_elem = {0., 0.,0.};
555 for (
int corn = 0; corn < 8; ++corn) {
556 for (
int c = 0; c < 3; ++c)
558 center_in_father_reference_elem[c] += corners_in_father_reference_elem_temp[corn][c]/8.;
562 double volume_in_father_reference_elem = double(1)/(cells_per_dim[0]*cells_per_dim[1]*cells_per_dim[2]);
565 in_father_reference_elem_corners, in_father_reference_elem_corner_indices.data());
568 OPM_THROW(std::logic_error,
"Entity has no father.");
577 auto ancestor = this->
father();
578 while (ancestor.hasFather()){
579 ancestor = ancestor.father();
581 assert(ancestor.level() == 0);
584 else if (!(pgrid_ -> leaf_to_level_cells_.empty())) {
587 const int&
level = pgrid_->leaf_to_level_cells_[this->
index()][0];
589 const int& levelElemIdx = pgrid_->leaf_to_level_cells_[this->
index()][1];
603 if (!(pgrid_ -> leaf_to_level_cells_.empty()))
605 const int& entityLevelIdx = pgrid_->leaf_to_level_cells_[this->
index()][1];
616 const auto& level_data = (*(pgrid_ -> level_data_ptr_))[
level()].get();
617 return level_data -> global_cell_[
getLevelElem().index()];
[ provides Dune::Grid ]
Definition CpGrid.hpp:203
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:118
bool operator==(const EntityRep &other) const
Equality operator.
Definition EntityRep.hpp:179
int index() const
The (positive) index of an entity.
Definition EntityRep.hpp:126
EntityRep()
Default constructor.
Definition EntityRep.hpp:104
Base class for EntityVariable and SignedEntityVariable.
Definition EntityRep.hpp:219
bool mightVanish() const
Indicates whether the entity may be removed in the next call to adapt().
Definition Entity.hpp:494
Entity< 0 > father() const
ONLY FOR CELLS (Entity<0>).
Definition Entity.hpp:511
unsigned int subEntities(const unsigned int cc) const
Return the number of all subentities of the entity of a given codimension cc.
Definition Entity.hpp:380
LeafIntersectionIterator ileafend() const
End leaf-iterator for the cell-cell intersections of this entity.
Definition Entity.hpp:342
EntitySeed seed() const
Return an entity seed (light-weight entity).
Definition Entity.hpp:159
int getLevelCartesianIdx() const
Get Cartesian Index in the level grid view where the Entity was born.
Definition Entity.hpp:614
HierarchicIterator hend(int) const
Iterator end over the children/beyond last child iterator.
Definition Entity.hpp:358
const Geometry & geometry() const
Return the geometry of the entity (does not depend on its orientation).
Definition Entity.hpp:394
Entity(const CpGridData &grid, int index_arg, bool orientation_arg)
Constructor taking a grid, entity index, and orientation.
Definition Entity.hpp:133
bool hasFather() const
ONLY FOR CELLS (Entity<0>) Check if the entity comes from an LGR, i.e., it has been created via refin...
Definition Entity.hpp:500
bool isValid() const
isValid method for EntitySeed
Definition Entity.hpp:434
Entity< 0 > getOrigin() const
Returns (1) oldest ancestor, i.e., oldest parent entity in the level-grid 0, if the entity was born i...
Definition Entity.hpp:573
bool operator==(const Entity &other) const
Equality.
Definition Entity.hpp:145
bool isRegular() const
Refinement is not defined for CpGrid.
Definition Entity.hpp:178
static constexpr int codimension
Definition Entity.hpp:86
Entity()
Constructor taking a grid and an integer entity representation.
Definition Entity.hpp:121
HierarchicIterator hbegin(int) const
Iterator begin over the children. [If requested, also over descendants more than one generation away....
Definition Entity.hpp:350
const Entity & impl() const
Access the actual implementation class behind Entity interface class.
Definition Entity.hpp:266
PartitionType partitionType() const
In serial run, the only partitionType() is InteriorEntity.
Definition Entity.hpp:365
Dune::cpgrid::Geometry< 3, 3 > geometryInFather() const
Return LocalGeometry representing the embedding of the entity into its father (when hasFather() is tr...
Definition Entity.hpp:531
Entity(int index_arg, bool orientation_arg)
Constructor taking a entity index, and orientation.
Definition Entity.hpp:139
LeafIntersectionIterator ileafbegin() const
Start leaf-iterator for the cell-cell intersections of this entity.
Definition Entity.hpp:335
Entity(const CpGridData &grid, EntityRep< codim > entityrep)
Constructor taking a grid and an entity representation.
Definition Entity.hpp:127
LevelIntersectionIterator ilevelend() const
End level-iterator for the cell-cell intersections of this entity.
Definition Entity.hpp:328
bool isNew() const
Returns true, if the entity has been created during the last call to adapt().
Definition Entity.hpp:477
Entity< 0 > getLevelElem() const
Get equivalent element on the level grid where the entity was born, if grid = leaf-grid-view....
Definition Entity.hpp:598
LevelIntersectionIterator ilevelbegin() const
Start level-iterator for the cell-cell intersections of this entity.
Definition Entity.hpp:321
GeometryType type() const
Return marker object (GeometryType object) representing the reference element of the entity.
Definition Entity.hpp:192
int level() const
Return the level of the entity in the grid hierarchy. Level = 0 represents the coarsest grid.
Definition Entity.hpp:442
bool isLeaf() const
Check if the entity is in the leafview.
Definition Entity.hpp:466
bool operator!=(const Entity &other) const
Inequality.
Definition Entity.hpp:151
bool hasBoundaryIntersections() const
Returns true if any of my intersections are on the boundary.
Definition Entity.hpp:421
Codim< cc >::Entity subEntity(int i) const
Obtain subentity.
This class encapsulates geometry for vertices, intersections, and cells.
Definition Geometry.hpp:76
Only needs to provide interface for doing nothing.
Definition Iterators.hpp:118
Definition Intersection.hpp:279
Iterator intended to be used as LeafIterator and LevelIterator (no difference due to no adaptivity) f...
Definition Iterators.hpp:60
Definition Indexsets.hpp:371
Copyright 2019 Equinor AS.
Definition GridPartitioning.cpp:673
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
Export supported entity types.
Definition Entity.hpp:97