Entity.hpp
Go to the documentation of this file.
1//===========================================================================
2//
3// File: Entity.hpp
4//
5// Created: Fri May 29 20:26:48 2009
6//
7// Author(s): Atgeirr F Rasmussen <atgeirr@sintef.no>
8// Brd Skaflestad <bard.skaflestad@sintef.no>
9// Antonella Ritorto <antonella.ritorto@opm-op.com>
10//
11// $Date$
12//
13// $Revision$
14//
15//===========================================================================
16
17/*
18 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
19 Copyright 2009, 2010, 2022 Equinor ASA.
20
21 This file is part of The Open Porous Media project (OPM).
22
23 OPM is free software: you can redistribute it and/or modify
24 it under the terms of the GNU General Public License as published by
25 the Free Software Foundation, either version 3 of the License, or
26 (at your option) any later version.
27
28 OPM is distributed in the hope that it will be useful,
29 but WITHOUT ANY WARRANTY; without even the implied warranty of
30 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
31 GNU General Public License for more details.
32
33 You should have received a copy of the GNU General Public License
34 along with OPM. If not, see <http://www.gnu.org/licenses/>.
35*/
36
37#ifndef OPM_ENTITY_HEADER
38#define OPM_ENTITY_HEADER
39
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>
44
47
48// To be able to test local and global ids of vertices
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>&);
54
55namespace Dune
56{
57namespace cpgrid
58{
59
60template <int codim>
61class Entity;
62
63namespace Impl {
64// Forward declaration of traits class for codimensions.
65template <int codim>
67
68// Specialization of traits class for codimension 0 (cells).
69template<>
70struct CodimTraits<0>
71{
73};
74
75#if DUNE_VERSION_LT(DUNE_GRID, 2, 11)
76// Specialization of traits class for codimension 1 (faces).
77// This is only needed for DUNE < 2.11, as the face entity type was needed for the VTK writer, which was fixed in DUNE 2.11.
78// See: https://gitlab.dune-project.org/core/dune-grid/-/merge_requests/793
79template<>
80struct CodimTraits<1>
81{
83};
84#endif
85
86// Specialization of traits class for codimension 3 (vertices).
87template<>
88struct CodimTraits<3>
89{
91};
92
93}
94
95template<int,int> class Geometry;
96template<int,PartitionIteratorType> class Iterator;
99class CpGridData;
100class LevelGlobalIdSet;
101
105template <int codim>
106class Entity : public EntityRep<codim>
107{
108 friend class LevelGlobalIdSet;
109 friend class GlobalIdSet;
110 friend class HierarchicIterator;
111 friend class CpGridData;
113 const std::vector<std::array<int,3>>&,
114 const std::vector<std::array<int,3>>&,
115 const std::vector<std::array<int,3>>&,
116 const std::vector<std::string>&);
117
118public:
121 static constexpr int codimension = codim;
122 static constexpr int dimension = 3;
123 static constexpr int mydimension = dimension - codimension;
124 static constexpr int dimensionworld = 3;
125
126 // the official DUNE names
128
130 template <int cd>
131 using Codim = typename Impl::CodimTraits<cd>;
132
133 typedef cpgrid::Geometry<3-codim,3> Geometry;
135
139
140 typedef double ctype;
141
146 // Entity(const CpGridData& grid, int entityrep)
147 // : EntityRep<codim>(entityrep), pgrid_(&grid)
148 // {
149 // }
150
153 : EntityRep<codim>(), pgrid_( 0 )
154 {
155 }
156
158 Entity(const CpGridData& grid, EntityRep<codim> entityrep)
159 : EntityRep<codim>(entityrep), pgrid_(&grid)
160 {
161 }
162
164 Entity(const CpGridData& grid, int index_arg, bool orientation_arg)
165 : EntityRep<codim>(index_arg, orientation_arg), pgrid_(&grid)
166 {
167 }
168
170 Entity(int index_arg, bool orientation_arg)
171 : EntityRep<codim>(index_arg, orientation_arg), pgrid_()
172 {
173 }
174
176 bool operator==(const Entity& other) const
177 {
178 return EntityRep<codim>::operator==(other) && pgrid_ == other.pgrid_;
179 }
180
182 bool operator!=(const Entity& other) const
183 {
184 return !operator==(other);
185 }
186
191 {
192 return EntitySeed( impl() );
193 }
194
196 const Geometry& geometry() const;
197
199 int level() const;
200
206 bool isLeaf() const;
207
209 bool isRegular() const
210 {
211 return true;
212 }
213
219 PartitionType partitionType() const;
220
223 GeometryType type() const
224 {
225 return Dune::GeometryTypes::cube(3 - codim);
226 }
227
229 unsigned int subEntities ( const unsigned int cc ) const;
230
233 template <int cc>
234 typename Codim<cc>::Entity subEntity(int i) const;
235
238
241
244
247
248
251
254
261 bool isNew() const;
262
268 bool mightVanish() const;
269
275 bool hasFather() const;
276
281
288
294
295 // Mimic Dune entity wrapper
297 const Entity& impl() const
298 {
299 return *this;
300 }
301
303 {
304 return *this;
305 }
306
309 bool isValid () const;
310
319
322
325
327
328protected:
330};
331
332// Reference element for a given entity
333template <int codim>
335{
336 return Dune::referenceElement<double,3-codim>(Dune::GeometryTypes::cube(3-codim));
337}
338
339
340} // namespace cpgrid
341} // namespace Dune
342
343// now we include the Iterators.hh We need to do this here because for hbegin/hend the compiler
344// needs to know the size of hierarchicIterator
345#include "Iterators.hpp"
346#include "Intersection.hpp"
347namespace Dune
348{
349namespace cpgrid
350{
351template<int codim>
353{
354 static_assert(codim == 0, "");
355 return LevelIntersectionIterator(*pgrid_, *this, false);
356}
357
358template<int codim>
360{
361 static_assert(codim == 0, "");
362 return LevelIntersectionIterator(*pgrid_, *this, true);
363}
364
365template<int codim>
367{
368 static_assert(codim == 0, "");
369 return LeafIntersectionIterator(*pgrid_, *this, false);
370}
371
372template<int codim>
374{
375 static_assert(codim == 0, "");
376 return LeafIntersectionIterator(*pgrid_, *this, true);
377}
378
379
380template<int codim>
382{
383 // Creates iterator with first child as target if there is one. Otherwise empty stack and target.
384 return HierarchicIterator(*this, maxLevel);
385}
386
388template<int codim>
390{
391 // Creates iterator with empty stack and target.
392 return HierarchicIterator(maxLevel);
393}
394
395template <int codim>
396PartitionType Entity<codim>::partitionType() const
397{
398 return pgrid_->partition_type_indicator_->getPartitionType(*this);
399}
400
401} // namespace cpgrid
402} // namespace Dune
403
404
406
407namespace Dune {
408namespace cpgrid {
409
410template<int codim>
411unsigned int Entity<codim>::subEntities ( const unsigned int cc ) const
412{
413 if (cc == codim) {
414 return 1;
415 }
416 else if ( codim == 0 ){ // Cell/element/Entity<0>
417 if ( cc == 3 ) { // Get number of corners of the element.
418 return 8;
419 }
420 }
421 DUNE_THROW(NotImplemented, "subEntities not implemented for this codimension");
422}
423
424template <int codim>
426{
427 return pgrid_->geomVector<codim>()[*this];
428}
429
430template <int codim>
431template <int cc>
433{
434 static_assert(codim == 0, "");
435 if (cc == 0) { // Cell/element/Entity<0>
436 assert(i == 0);
437 typename Codim<cc>::Entity se(*pgrid_, this->index(), this->orientation());
438 return se;
439 } else if (cc == 3) { // Corner/Entity<3>
440 assert(i >= 0 && i < 8);
441 int corner_index = pgrid_->cell_to_point_[this->index()][i];
442 typename Codim<cc>::Entity se(*pgrid_, corner_index, true);
443 return se;
444 }
445 else {
446 OPM_THROW(std::runtime_error,
447 "No subentity exists of codimension " + std::to_string(cc));
448 }
449}
450
451template <int codim>
453{
454 // Copied implementation from EntityDefaultImplementation,
455 // except for not checking LevelIntersectionIterators.
456 typedef LeafIntersectionIterator Iter;
457 Iter end = ileafend();
458 for (Iter it = ileafbegin(); it != end; ++it) {
459 if (it->boundary()) return true;
460 }
461 return false;
462}
463
464template <int codim>
466{
467 return pgrid_ ? EntityRep<codim>::index() < pgrid_->size(codim) : false;
468}
469
470
471// level() It simply returns the level of the entity in the grid hierarchy.
472template <int codim>
474{
475 // Parallel and LGRs:
476 // If the grid has been distributed and - after that - LGRs have been added, then level_data_ptr_
477 // points at distributed_data_ which has size > 1.
478 //
479 // Serial and LGRs:
480 // If the grid is not distributed and there LGRs have been added, level_data_ptr_ points at data_ which
481 // has size > 1.
482 //
483 // - leaf_to_level_cells_ is non-empty only on the leaf grid view of a grid that has been refined.
484 // - level_ is set equal to zero when instantiating a pgrid, and rewriten when it corresponds to a refined level grid.
485 return pgrid_->leaf_to_level_cells_.empty()? pgrid_->level_ : pgrid_->leaf_to_level_cells_[this-> index()][0];
486}
487
488// isLeaf()
489// - if distributed_data_ is empty: an element is a leaf <-> hbegin and hend return the same iterator. Then,
490// *cells from level 0 (coarse grid) that are not parents, are Leaf.
491// *cells from any LGR are Leaf, since they do not children (nested refinement not supported).
492// - if distrubuted_data_ is NOT empty: there may be children on a different process. Therefore,
493// isLeaf() returns true <-> the element is a leaf entity of the global refinement hierarchy. Equivalently,
494// it can be checked whether parent_to_children_cells_ is empty.
495
496template<int codim>
498{
499 if (pgrid_ -> parent_to_children_cells_.empty()){ // Grid cells without nested refinement
500 return true;
501 }
502 else {
503 return (std::get<0>((pgrid_ -> parent_to_children_cells_)[this-> index()]) == -1); // Not involved in any LGR
504 }
505}
506
507template<int codim>
509{
510 int data_count = pgrid_->level_data_ptr_->size();
511 if (data_count == 1) { // CpGrid has only level zero grid.
512 return false; // no element is new
513 }
514 else {
515 assert(data_count >= 2); // Level zero and leaf grids, plus at least one refined level grid
516 // In nested refinement, an element created during the last refinement step
517 // may itself be refined further within the same step at a higher level.
518 // To handle this case, we also check whether the element is a leaf
519 // (i.e., it has no children).
520 return isLeaf() && (pgrid_->refinement_max_level_ == data_count - 2); // minus "level zero" and "leaf"
521 }
522}
523
524template<int codim>
526{
527 return false;
528}
529
530template<int codim>
532{
533 if ((pgrid_ -> child_to_parent_cells_.empty()) || (pgrid_ -> child_to_parent_cells_[this->index()][0] == -1)){
534 return false;
535 }
536 else{
537 return true;
538 }
539}
540
541template<int codim>
543{
544 if (this->hasFather()){
545 const int& coarser_level = pgrid_ -> child_to_parent_cells_[this->index()][0];
546 const int& parent_cell_index = pgrid_ -> child_to_parent_cells_[this->index()][1];
547 return Entity<0>( *((*(pgrid_ -> level_data_ptr_))[coarser_level].get()), parent_cell_index, true);
548 }
549 else{
550 OPM_THROW(std::logic_error, "Entity has no father.");
551 }
552}
553
554template<int codim>
556{
557 return pgrid_ -> cell_to_idxInParentCell_[this->index()];
558}
559
560
561template<int codim>
563{
564 if (!(this->hasFather())){
565 OPM_THROW(std::logic_error, "Entity has no father.");
566 }
567
568 // Indices of corners in entity's geometry in father reference element.
569 static constexpr std::array<int,8> in_father_reference_elem_corner_indices = {0,1,2,3,4,5,6,7};
570 // 'static': The returned object Geometry<3,3> stores a pointer to in_father_reference_elem_corner_indices. Therefore,
571 // this variable is declared static to prolongate its lifetime beyond this function (static storage duration).
572
573 auto idx_in_parent_cell = pgrid_ -> cell_to_idxInParentCell_[this->index()];
574 if (idx_in_parent_cell !=-1) {
575 const auto& cells_per_dim = (*(pgrid_ -> level_data_ptr_))[this->level()] -> cells_per_dim_;
576 const auto& auxArr = pgrid_ -> getReferenceRefinedCorners(idx_in_parent_cell, cells_per_dim);
577 FieldVector<double, 3> corners_in_father_reference_elem_temp[8] =
578 { auxArr[0], auxArr[1], auxArr[2], auxArr[3], auxArr[4], auxArr[5], auxArr[6], auxArr[7]};
579 auto in_father_reference_elem_corners = std::make_shared<EntityVariable<cpgrid::Geometry<0, 3>, 3>>();
580 EntityVariableBase<cpgrid::Geometry<0, 3>>& mutable_in_father_reference_elem_corners = *in_father_reference_elem_corners;
581 // Assign the corners. Make use of the fact that pointers behave like iterators.
582 mutable_in_father_reference_elem_corners.assign(corners_in_father_reference_elem_temp,
583 corners_in_father_reference_elem_temp + 8);
584 // Compute the center of the 'local-entity'.
585 Dune::FieldVector<double, 3> center_in_father_reference_elem = {0., 0.,0.};
586 for (int corn = 0; corn < 8; ++corn) {
587 for (int c = 0; c < 3; ++c)
588 {
589 center_in_father_reference_elem[c] += corners_in_father_reference_elem_temp[corn][c]/8.;
590 }
591 }
592 // Compute the volume of the 'local-entity'.
593 double volume_in_father_reference_elem = double(1)/(cells_per_dim[0]*cells_per_dim[1]*cells_per_dim[2]);
594 // Construct (and return) the Geometry<3,3> of 'child-cell in the reference element of its father (unit cube)'.
595 return Dune::cpgrid::Geometry<3,3>(center_in_father_reference_elem, volume_in_father_reference_elem,
596 in_father_reference_elem_corners, in_father_reference_elem_corner_indices.data());
597 }
598 else {
599 OPM_THROW(std::logic_error, "Entity has no father.");
600 }
601}
602
603template<int codim>
605{
606 if (hasFather()) { // Entit is a refined cell belonging to
607 // a refined level grid (level>0) or to the leaf grid.
608 auto ancestor = this->father();
609 while (ancestor.hasFather()){
610 ancestor = ancestor.father();
611 }
612 assert(ancestor.level() == 0);
613 return ancestor;
614 }
615 else if (!(pgrid_ -> leaf_to_level_cells_.empty())) { // Entity is a coarse cell
616 // (born in level zero grid, never involved in refinement) belonging to the leaf grid.
617 // leaf_to_level_cells_ [leaf idx] = { level where entity was born, cell idx in that level}
618 const int& level = pgrid_->leaf_to_level_cells_[this->index()][0];
619 assert(level == 0);
620 const int& levelElemIdx = pgrid_->leaf_to_level_cells_[this->index()][1];
621 return Dune::cpgrid::Entity<0>( *((*(pgrid_ -> level_data_ptr_))[level].get()), levelElemIdx, true);
622 }
623 else { // Entity belongs to level zero grid.
624 return *this;
625 }
626}
627
628template<int codim>
630{
631 // Check that the element belongs to the leaf grid view
632 // This is needed to get the index of the element in the level it was born.
633 // leaf_to_level_cells_ [leaf idx] = {level where the entity was born, equivalent cell idx in that level}
634 if (!(pgrid_ -> leaf_to_level_cells_.empty())) // entity on the LeafGridView
635 {
636 const int& entityLevelIdx = pgrid_->leaf_to_level_cells_[this->index()][1];
637 return Dune::cpgrid::Entity<0>( *((*(pgrid_ -> level_data_ptr_))[this->level()].get()), entityLevelIdx, true);
638 }
639 else {
640 return *this;
641 }
642}
643
644template<int codim>
646{
647 const auto& level_data = (*(pgrid_ -> level_data_ptr_))[level()].get();
648 return level_data -> global_cell_[getLevelElem().index()];
649}
650
651} // namespace cpgrid
652} // namespace Dune
653
654
655#endif // OPM_ENTITY_HEADER
void refinePatch_and_check(Dune::CpGrid &, const std::vector< std::array< int, 3 > > &, const std::vector< std::array< int, 3 > > &, const std::vector< std::array< int, 3 > > &, const std::vector< std::string > &)
[ provides Dune::Grid ]
Definition: CpGrid.hpp:203
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:118
bool mightVanish() const
Indicates whether the entity may be removed in the next call to adapt().
Definition: Entity.hpp:525
Entity< 0 > father() const
ONLY FOR CELLS (Entity<0>). Get the father Entity (in the level-grid the father cell was born),...
Definition: Entity.hpp:542
const CpGridData * pgrid_
Definition: Entity.hpp:329
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:411
LeafIntersectionIterator ileafend() const
End leaf-iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:373
EntitySeed seed() const
Return an entity seed (light-weight entity). EntitySeed objects are used to obtain an Entity back whe...
Definition: Entity.hpp:190
int getLevelCartesianIdx() const
Get Cartesian Index in the level grid view where the Entity was born.
Definition: Entity.hpp:645
HierarchicIterator hend(int) const
Iterator end over the children/beyond last child iterator.
Definition: Entity.hpp:389
const Geometry & geometry() const
Return the geometry of the entity (does not depend on its orientation).
Definition: Entity.hpp:425
cpgrid::IntersectionIterator LeafIntersectionIterator
Definition: Entity.hpp:136
double ctype
Definition: Entity.hpp:140
Entity(const CpGridData &grid, int index_arg, bool orientation_arg)
Constructor taking a grid, entity index, and orientation.
Definition: Entity.hpp:164
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:531
bool isValid() const
Definition: Entity.hpp:465
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:604
Geometry LocalGeometry
Definition: Entity.hpp:134
cpgrid::HierarchicIterator HierarchicIterator
Definition: Entity.hpp:138
bool operator==(const Entity &other) const
Equality.
Definition: Entity.hpp:176
bool isRegular() const
Refinement is not defined for CpGrid.
Definition: Entity.hpp:209
static constexpr int codimension
Definition: Entity.hpp:121
static constexpr int mydimension
Definition: Entity.hpp:123
int getIdxInParentCell() const
Definition: Entity.hpp:555
Entity()
Constructor creating empty entity.
Definition: Entity.hpp:152
HierarchicIterator hbegin(int) const
Iterator begin over the children. [If requested, also over descendants more than one generation away....
Definition: Entity.hpp:381
cpgrid::Geometry< 3-codim, 3 > Geometry
Definition: Entity.hpp:133
cpgrid::IntersectionIterator LevelIntersectionIterator
Definition: Entity.hpp:137
static constexpr int dimension
Definition: Entity.hpp:122
const Entity & impl() const
Access the actual implementation class behind Entity interface class.
Definition: Entity.hpp:297
PartitionType partitionType() const
In serial run, the only partitionType() is InteriorEntity. Only needed when distributed_data_ is not ...
Definition: Entity.hpp:396
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:562
typename Impl::CodimTraits< cd > Codim
Export supported entity types.
Definition: Entity.hpp:131
Entity(int index_arg, bool orientation_arg)
Constructor taking a entity index, and orientation.
Definition: Entity.hpp:170
LeafIntersectionIterator ileafbegin() const
Start leaf-iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:366
static constexpr int dimensionworld
Definition: Entity.hpp:124
Entity(const CpGridData &grid, EntityRep< codim > entityrep)
Constructor taking a grid and an entity representation.
Definition: Entity.hpp:158
Entity & impl()
Definition: Entity.hpp:302
LevelIntersectionIterator ilevelend() const
End level-iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:359
bool isNew() const
Returns true, if the entity has been created during the last call to adapt().
Definition: Entity.hpp:508
Entity< 0 > getLevelElem() const
Get equivalent element on the level grid where the entity was born, if grid = leaf-grid-view....
Definition: Entity.hpp:629
LevelIntersectionIterator ilevelbegin() const
Start level-iterator for the cell-cell intersections of this entity.
Definition: Entity.hpp:352
GeometryType type() const
Return marker object (GeometryType object) representing the reference element of the entity....
Definition: Entity.hpp:223
int level() const
Return the level of the entity in the grid hierarchy. Level = 0 represents the coarsest grid.
Definition: Entity.hpp:473
bool isLeaf() const
Check if the entity is in the leafview.
Definition: Entity.hpp:497
bool operator!=(const Entity &other) const
Inequality.
Definition: Entity.hpp:182
Entity EntitySeed
Definition: Entity.hpp:127
bool hasBoundaryIntersections() const
Definition: Entity.hpp:452
Codim< cc >::Entity subEntity(int i) const
Obtain subentity. Example: If cc = 3 and i = 5, it returns the 5th corner/vertex of the entity.
Represents an entity of a given codim, with positive or negative orientation.
Definition: EntityRep.hpp:99
bool operator==(const EntityRep &other) const
Equality operator.
Definition: EntityRep.hpp:179
int index() const
The (positive) index of an entity. Not a Dune interface method.
Definition: EntityRep.hpp:126
Base class for EntityVariable and SignedEntityVariable. Forwards a restricted subset of the std::vect...
Definition: EntityRep.hpp:219
Definition: Geometry.hpp:76
The global id set for Dune.
Definition: Indexsets.hpp:483
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Intersection.hpp:279
Definition: Iterators.hpp:60
Definition: Indexsets.hpp:367
auto referenceElement(const Entity< codim > &)
Definition: Entity.hpp:334
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
Definition: Entity.hpp:66