dune-mmesh (1.4)

intersectioniterator.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_MMESH_INTERFACE_INTERSECTIONITERATOR_HH
4#define DUNE_MMESH_INTERFACE_INTERSECTIONITERATOR_HH
5
10// Dune includes
11#include <dune/grid/common/intersection.hh>
12
13// Dune MMesh includes
16
17namespace Dune
18{
19
28 template<class GridImp>
30 {
31
32 enum {dimension=GridImp::dimension};
33
34 enum {dimensionworld=GridImp::dimensionworld};
35
36 // The type used to store coordinates
37 typedef typename GridImp::ctype ctype;
38
39 typedef typename GridImp::template MMeshInterfaceEntity<0> MMeshInterfaceEntity;
40
41 public:
42 typedef Dune::Intersection<const GridImp, Dune::MMeshInterfaceGridLeafIntersection<GridImp> > Intersection;
43
46
49 const MMeshInterfaceEntity& hostEntity)
50 : grid_(grid)
51 , hostEntity_(hostEntity)
52 , i_( 0 )
53 , nbIdx_( 0 )
54 {
55 const auto& indexSet = grid_->leafIndexSet();
56 const auto& cgalIndex = MMeshInterfaceImpl::computeCGALIndices<MMeshInterfaceEntity, dimension>( hostEntity );
57
58 for( int d = 0; d < dimension+1; ++d )
59 {
60 std::array< std::size_t, dimension > ids;
61 if constexpr (dimension == 1)
62 {
63 ids[0] = indexSet.vertexIndexMap().at( hostEntity.first->vertex( cgalIndex[d] )->info().id );
64 }
65 else // dim == 2
66 {
67 ids[0] = indexSet.vertexIndexMap().at( hostEntity.first->vertex( cgalIndex[d==2 ? 1 : 0] )->info().id );
68 ids[1] = indexSet.vertexIndexMap().at( hostEntity.first->vertex( cgalIndex[d==0 ? 1 : 2] )->info().id );
69 }
70
71 try {
72 std::sort(ids.begin(), ids.end());
73 maxNbIdx_[d] = (int)indexSet.indexMap().at( ids ).size() - 2;
74 } catch (std::exception &e) {
75 DUNE_THROW(InvalidStateException, e.what());
76 }
77 }
78
79 while( skip() )
80 increment();
81 }
82
85 const MMeshInterfaceEntity& hostEntity,
86 bool endDummy)
87 : grid_(grid)
88 , hostEntity_(hostEntity)
89 , i_( dimension + 1 )
90 , nbIdx_( 0 )
91 {}
92
95 return hostEntity_ == other.hostEntity_
96 && i_ == other.i_
97 && nbIdx_ == other.nbIdx_;
98 }
99
102 {
103 if ( maxNbIdx_[i_] != -1 and nbIdx_ < maxNbIdx_[i_] )
104 nbIdx_++;
105 else
106 {
107 ++i_;
108 nbIdx_ = 0;
109 }
110
111 if (skip())
112 increment();
113 }
114
116 Intersection dereference() const {
117 return MMeshInterfaceGridLeafIntersection<GridImp>(grid_, hostEntity_, i_, nbIdx_);
118 }
119
120 private:
121 bool skip()
122 {
123 if (i_ == dimension + 1)
124 return false;
125 if (grid_->entity(hostEntity_).partitionType() == GhostEntity)
126 return maxNbIdx_[i_] == -1 || dereference().outside().partitionType() == GhostEntity;
127 return false;
128 }
129
130 const GridImp* grid_;
131 MMeshInterfaceEntity hostEntity_;
132 std::size_t i_;
133 std::size_t nbIdx_;
134 std::array<int, dimension+1> maxNbIdx_;
135 };
136
137} // namespace Dune
138
139#endif
Iterator over all element neighborsMesh entities of codimension 0 ("elements") allow to visit all nei...
Definition: intersectioniterator.hh:30
MMeshInterfaceGridLeafIntersectionIterator(const GridImp *grid, const MMeshInterfaceEntity &hostEntity)
constructor for (begin) iterator
Definition: intersectioniterator.hh:48
MMeshInterfaceGridLeafIntersectionIterator(const GridImp *grid, const MMeshInterfaceEntity &hostEntity, bool endDummy)
constructor for end iterator
Definition: intersectioniterator.hh:84
Intersection dereference() const
dereferencing
Definition: intersectioniterator.hh:116
void increment()
prefix increment
Definition: intersectioniterator.hh:101
bool equals(const MMeshInterfaceGridLeafIntersectionIterator &other) const
returns if iterators reference same intersection
Definition: intersectioniterator.hh:94
MMeshInterfaceGridLeafIntersectionIterator()
default constructor
Definition: intersectioniterator.hh:45
An intersection with a leaf neighbor elementMesh entities of codimension 0 ("elements") allow to visi...
Definition: intersections.hh:37
The MMeshInterfaceGridEntity class.
The MMeshInterfaceGridLeafIntersection and MMeshLevelIntersection classes.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 7, 22:57, 2025)