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,int> class Geometry;
61template<int,PartitionIteratorType> class Iterator;
62class IntersectionIterator;
63class HierarchicIterator;
64class CpGridData;
65class LevelGlobalIdSet;
66
70template <int codim>
71class Entity : public EntityRep<codim>
72{
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>&);
82
83public:
86 static constexpr int codimension = codim;
87 static constexpr int dimension = 3;
88 static constexpr int mydimension = dimension - codimension;
89 static constexpr int dimensionworld = 3;
90
91 // the official DUNE names
93
95 template <int cd>
96 struct Codim
97 {
99 };
100
101
102 typedef cpgrid::Geometry<3-codim,3> Geometry;
104
108
109 typedef double ctype;
110
115 // Entity(const CpGridData& grid, int entityrep)
116 // : EntityRep<codim>(entityrep), pgrid_(&grid)
117 // {
118 // }
119
122 : EntityRep<codim>(), pgrid_( 0 )
123 {
124 }
125
127 Entity(const CpGridData& grid, EntityRep<codim> entityrep)
128 : EntityRep<codim>(entityrep), pgrid_(&grid)
129 {
130 }
131
133 Entity(const CpGridData& grid, int index_arg, bool orientation_arg)
134 : EntityRep<codim>(index_arg, orientation_arg), pgrid_(&grid)
135 {
136 }
137
139 Entity(int index_arg, bool orientation_arg)
140 : EntityRep<codim>(index_arg, orientation_arg), pgrid_()
141 {
142 }
143
145 bool operator==(const Entity& other) const
146 {
147 return EntityRep<codim>::operator==(other) && pgrid_ == other.pgrid_;
148 }
149
151 bool operator!=(const Entity& other) const
152 {
153 return !operator==(other);
154 }
155
160 {
161 return EntitySeed( impl() );
162 }
163
165 const Geometry& geometry() const;
166
168 int level() const;
169
175 bool isLeaf() const;
176
178 bool isRegular() const
179 {
180 return true;
181 }
182
188 PartitionType partitionType() const;
189
192 GeometryType type() const
193 {
194 return Dune::GeometryTypes::cube(3 - codim);
195 }
196
198 unsigned int subEntities ( const unsigned int cc ) const;
199
202 template <int cc>
203 typename Codim<cc>::Entity subEntity(int i) const;
204
207
210
213
216
217
220
223
230 bool isNew() const;
231
237 bool mightVanish() const;
238
244 bool hasFather() const;
245
250
257
263
264 // Mimic Dune entity wrapper
266 const Entity& impl() const
267 {
268 return *this;
269 }
270
272 {
273 return *this;
274 }
275
278 bool isValid () const;
279
288
291
294
296
297protected:
299};
300
301// Reference element for a given entity
302template <int codim>
304{
305 return Dune::referenceElement<double,3-codim>(Dune::GeometryTypes::cube(3-codim));
306}
307
308
309} // namespace cpgrid
310} // namespace Dune
311
312// now we include the Iterators.hh We need to do this here because for hbegin/hend the compiler
313// needs to know the size of hierarchicIterator
314#include "Iterators.hpp"
315#include "Intersection.hpp"
316namespace Dune
317{
318namespace cpgrid
319{
320template<int codim>
322{
323 static_assert(codim == 0, "");
324 return LevelIntersectionIterator(*pgrid_, *this, false);
325}
326
327template<int codim>
329{
330 static_assert(codim == 0, "");
331 return LevelIntersectionIterator(*pgrid_, *this, true);
332}
333
334template<int codim>
336{
337 static_assert(codim == 0, "");
338 return LeafIntersectionIterator(*pgrid_, *this, false);
339}
340
341template<int codim>
343{
344 static_assert(codim == 0, "");
345 return LeafIntersectionIterator(*pgrid_, *this, true);
346}
347
348
349template<int codim>
351{
352 // Creates iterator with first child as target if there is one. Otherwise empty stack and target.
353 return HierarchicIterator(*this, maxLevel);
354}
355
357template<int codim>
359{
360 // Creates iterator with empty stack and target.
361 return HierarchicIterator(maxLevel);
362}
363
364template <int codim>
365PartitionType Entity<codim>::partitionType() const
366{
367 return pgrid_->partition_type_indicator_->getPartitionType(*this);
368}
369
370} // namespace cpgrid
371} // namespace Dune
372
373
375
376namespace Dune {
377namespace cpgrid {
378
379template<int codim>
380unsigned int Entity<codim>::subEntities ( const unsigned int cc ) const
381{
382 if (cc == codim) {
383 return 1;
384 }
385 else if ( codim == 0 ){ // Cell/element/Entity<0>
386 if ( cc == 3 ) { // Get number of corners of the element.
387 return 8;
388 }
389 }
390 return 0;
391}
392
393template <int codim>
395{
396 return pgrid_->geomVector<codim>()[*this];
397}
398
399template <int codim>
400template <int cc>
402{
403 static_assert(codim == 0, "");
404 if (cc == 0) { // Cell/element/Entity<0>
405 assert(i == 0);
406 typename Codim<cc>::Entity se(*pgrid_, this->index(), this->orientation());
407 return se;
408 } else if (cc == 3) { // Corner/Entity<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);
412 return se;
413 }
414 else {
415 OPM_THROW(std::runtime_error,
416 "No subentity exists of codimension " + std::to_string(cc));
417 }
418}
419
420template <int codim>
422{
423 // Copied implementation from EntityDefaultImplementation,
424 // except for not checking LevelIntersectionIterators.
425 typedef LeafIntersectionIterator Iter;
426 Iter end = ileafend();
427 for (Iter it = ileafbegin(); it != end; ++it) {
428 if (it->boundary()) return true;
429 }
430 return false;
431}
432
433template <int codim>
435{
436 return pgrid_ ? EntityRep<codim>::index() < pgrid_->size(codim) : false;
437}
438
439
440// level() It simply returns the level of the entity in the grid hierarchy.
441template <int codim>
443{
444 // Parallel and LGRs:
445 // If the grid has been distributed and - after that - LGRs have been added, then level_data_ptr_
446 // points at distributed_data_ which has size > 1.
447 //
448 // Serial and LGRs:
449 // If the grid is not distributed and there LGRs have been added, level_data_ptr_ points at data_ which
450 // has size > 1.
451 //
452 // - leaf_to_level_cells_ is non-empty only on the leaf grid view of a grid that has been refined.
453 // - level_ is set equal to zero when instantiating a pgrid, and rewriten when it corresponds to a refined level grid.
454 return pgrid_->leaf_to_level_cells_.empty()? pgrid_->level_ : pgrid_->leaf_to_level_cells_[this-> index()][0];
455}
456
457// isLeaf()
458// - if distributed_data_ is empty: an element is a leaf <-> hbegin and hend return the same iterator. Then,
459// *cells from level 0 (coarse grid) that are not parents, are Leaf.
460// *cells from any LGR are Leaf, since they do not children (nested refinement not supported).
461// - if distrubuted_data_ is NOT empty: there may be children on a different process. Therefore,
462// isLeaf() returns true <-> the element is a leaf entity of the global refinement hierarchy. Equivalently,
463// it can be checked whether parent_to_children_cells_ is empty.
464
465template<int codim>
467{
468 if (pgrid_ -> parent_to_children_cells_.empty()){ // Grid cells without nested refinement
469 return true;
470 }
471 else {
472 return (std::get<0>((pgrid_ -> parent_to_children_cells_)[this-> index()]) == -1); // Not involved in any LGR
473 }
474}
475
476template<int codim>
478{
479 int data_count = pgrid_->level_data_ptr_->size();
480 if (data_count == 1) { // CpGrid has only level zero grid.
481 return false; // no element is new
482 }
483 else {
484 assert(data_count >= 2); // Level zero and leaf grids, plus at least one refined level grid
485 // In nested refinement, an element created during the last refinement step
486 // may itself be refined further within the same step at a higher level.
487 // To handle this case, we also check whether the element is a leaf
488 // (i.e., it has no children).
489 return isLeaf() && (pgrid_->refinement_max_level_ == data_count - 2); // minus "level zero" and "leaf"
490 }
491}
492
493template<int codim>
495{
496 return false;
497}
498
499template<int codim>
501{
502 if ((pgrid_ -> child_to_parent_cells_.empty()) || (pgrid_ -> child_to_parent_cells_[this->index()][0] == -1)){
503 return false;
504 }
505 else{
506 return true;
507 }
508}
509
510template<int codim>
512{
513 if (this->hasFather()){
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);
517 }
518 else{
519 OPM_THROW(std::logic_error, "Entity has no father.");
520 }
521}
522
523template<int codim>
525{
526 return pgrid_ -> cell_to_idxInParentCell_[this->index()];
527}
528
529
530template<int codim>
532{
533 if (!(this->hasFather())){
534 OPM_THROW(std::logic_error, "Entity has no father.");
535 }
536
537 // Indices of corners in entity's geometry in father reference element.
538 static constexpr std::array<int,8> in_father_reference_elem_corner_indices = {0,1,2,3,4,5,6,7};
539 // 'static': The returned object Geometry<3,3> stores a pointer to in_father_reference_elem_corner_indices. Therefore,
540 // this variable is declared static to prolongate its lifetime beyond this function (static storage duration).
541
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>>();
549 EntityVariableBase<cpgrid::Geometry<0, 3>>& mutable_in_father_reference_elem_corners = *in_father_reference_elem_corners;
550 // Assign the corners. Make use of the fact that pointers behave like iterators.
551 mutable_in_father_reference_elem_corners.assign(corners_in_father_reference_elem_temp,
552 corners_in_father_reference_elem_temp + 8);
553 // Compute the center of the 'local-entity'.
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)
557 {
558 center_in_father_reference_elem[c] += corners_in_father_reference_elem_temp[corn][c]/8.;
559 }
560 }
561 // Compute the volume of the 'local-entity'.
562 double volume_in_father_reference_elem = double(1)/(cells_per_dim[0]*cells_per_dim[1]*cells_per_dim[2]);
563 // Construct (and return) the Geometry<3,3> of 'child-cell in the reference element of its father (unit cube)'.
564 return Dune::cpgrid::Geometry<3,3>(center_in_father_reference_elem, volume_in_father_reference_elem,
565 in_father_reference_elem_corners, in_father_reference_elem_corner_indices.data());
566 }
567 else {
568 OPM_THROW(std::logic_error, "Entity has no father.");
569 }
570}
571
572template<int codim>
574{
575 if (hasFather()) { // Entit is a refined cell belonging to
576 // a refined level grid (level>0) or to the leaf grid.
577 auto ancestor = this->father();
578 while (ancestor.hasFather()){
579 ancestor = ancestor.father();
580 }
581 assert(ancestor.level() == 0);
582 return ancestor;
583 }
584 else if (!(pgrid_ -> leaf_to_level_cells_.empty())) { // Entity is a coarse cell
585 // (born in level zero grid, never involved in refinement) belonging to the leaf grid.
586 // leaf_to_level_cells_ [leaf idx] = { level where entity was born, cell idx in that level}
587 const int& level = pgrid_->leaf_to_level_cells_[this->index()][0];
588 assert(level == 0);
589 const int& levelElemIdx = pgrid_->leaf_to_level_cells_[this->index()][1];
590 return Dune::cpgrid::Entity<0>( *((*(pgrid_ -> level_data_ptr_))[level].get()), levelElemIdx, true);
591 }
592 else { // Entity belongs to level zero grid.
593 return *this;
594 }
595}
596
597template<int codim>
599{
600 // Check that the element belongs to the leaf grid view
601 // This is needed to get the index of the element in the level it was born.
602 // leaf_to_level_cells_ [leaf idx] = {level where the entity was born, equivalent cell idx in that level}
603 if (!(pgrid_ -> leaf_to_level_cells_.empty())) // entity on the LeafGridView
604 {
605 const int& entityLevelIdx = pgrid_->leaf_to_level_cells_[this->index()][1];
606 return Dune::cpgrid::Entity<0>( *((*(pgrid_ -> level_data_ptr_))[this->level()].get()), entityLevelIdx, true);
607 }
608 else {
609 return *this;
610 }
611}
612
613template<int codim>
615{
616 const auto& level_data = (*(pgrid_ -> level_data_ptr_))[level()].get();
617 return level_data -> global_cell_[getLevelElem().index()];
618}
619
620} // namespace cpgrid
621} // namespace Dune
622
623
624#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
Definition: Entity.hpp:72
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>). Get the father Entity (in the level-grid the father cell was born),...
Definition: Entity.hpp:511
const CpGridData * pgrid_
Definition: Entity.hpp:298
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). EntitySeed objects are used to obtain an Entity back whe...
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
cpgrid::IntersectionIterator LeafIntersectionIterator
Definition: Entity.hpp:105
double ctype
Definition: Entity.hpp:109
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
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
Geometry LocalGeometry
Definition: Entity.hpp:103
cpgrid::HierarchicIterator HierarchicIterator
Definition: Entity.hpp:107
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
static constexpr int mydimension
Definition: Entity.hpp:88
int getIdxInParentCell() const
Definition: Entity.hpp:524
Entity()
Constructor creating empty entity.
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
cpgrid::Geometry< 3-codim, 3 > Geometry
Definition: Entity.hpp:102
cpgrid::IntersectionIterator LevelIntersectionIterator
Definition: Entity.hpp:106
static constexpr int dimension
Definition: Entity.hpp:87
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. Only needed when distributed_data_ is not ...
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
static constexpr int dimensionworld
Definition: Entity.hpp:89
Entity(const CpGridData &grid, EntityRep< codim > entityrep)
Constructor taking a grid and an entity representation.
Definition: Entity.hpp:127
Entity & impl()
Definition: Entity.hpp:271
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
Entity EntitySeed
Definition: Entity.hpp:92
bool hasBoundaryIntersections() const
Definition: Entity.hpp:421
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:487
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:118
Definition: Intersection.hpp:279
Definition: Indexsets.hpp:371
auto referenceElement(const Entity< codim > &)
Definition: Entity.hpp:303
The namespace Dune is the main namespace for all Dune code.
Definition: common/CartesianIndexMapper.hpp:10
Export supported entity types.
Definition: Entity.hpp:97