dune-mmesh (unstable)

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
27template <class GridImp>
29 enum { dimension = GridImp::dimension };
30
31 enum { dimensionworld = GridImp::dimensionworld };
32
33 // The type used to store coordinates
34 typedef typename GridImp::ctype ctype;
35
36 typedef typename GridImp::template MMeshInterfaceEntity<0>
37 MMeshInterfaceEntity;
38
39 public:
40 typedef Dune::Intersection<const GridImp,
42 Intersection;
43
46
49 const GridImp* grid, const MMeshInterfaceEntity& hostEntity)
50 : grid_(grid), hostEntity_(hostEntity), i_(0), nbIdx_(0) {
51 const auto& indexSet = grid_->leafIndexSet();
52 const auto& cgalIndex =
53 MMeshInterfaceImpl::computeCGALIndices<MMeshInterfaceEntity, dimension>(
54 hostEntity);
55
56 for (int d = 0; d < dimension + 1; ++d) {
57 std::array<std::size_t, dimension> ids;
58 if constexpr (dimension == 1) {
59 ids[0] = indexSet.vertexIndexMap().at(
60 hostEntity.first->vertex(cgalIndex[d])->info().id);
61 } else // dim == 2
62 {
63 ids[0] = indexSet.vertexIndexMap().at(
64 hostEntity.first->vertex(cgalIndex[d == 2 ? 1 : 0])->info().id);
65 ids[1] = indexSet.vertexIndexMap().at(
66 hostEntity.first->vertex(cgalIndex[d == 0 ? 1 : 2])->info().id);
67 }
68
69 try {
70 std::sort(ids.begin(), ids.end());
71 maxNbIdx_[d] = (int)indexSet.indexMap().at(ids).size() - 2;
72 } catch (std::exception& e) {
73 DUNE_THROW(InvalidStateException, e.what());
74 }
75 }
76
77 while (skip()) increment();
78 }
79
82 const GridImp* grid, const MMeshInterfaceEntity& hostEntity,
83 bool endDummy)
84 : grid_(grid), hostEntity_(hostEntity), i_(dimension + 1), nbIdx_(0) {}
85
88 return hostEntity_ == other.hostEntity_ && i_ == other.i_ &&
89 nbIdx_ == other.nbIdx_;
90 }
91
93 void increment() {
94 if (maxNbIdx_[i_] != -1 and nbIdx_ < maxNbIdx_[i_])
95 nbIdx_++;
96 else {
97 ++i_;
98 nbIdx_ = 0;
99 }
100
101 if (skip()) increment();
102 }
103
105 Intersection dereference() const {
106 return MMeshInterfaceGridLeafIntersection<GridImp>(grid_, hostEntity_, i_,
107 nbIdx_);
108 }
109
110 private:
111 bool skip() {
112 if (i_ == dimension + 1) return false;
113 if (grid_->entity(hostEntity_).partitionType() == GhostEntity)
114 return maxNbIdx_[i_] == -1 ||
115 dereference().outside().partitionType() == GhostEntity;
116 return false;
117 }
118
119 const GridImp* grid_;
120 MMeshInterfaceEntity hostEntity_;
121 std::size_t i_;
122 std::size_t nbIdx_;
123 std::array<int, dimension + 1> maxNbIdx_;
124};
125
126} // namespace Dune
127
128#endif
Iterator over all element neighborsMesh entities of codimension 0 ("elements") allow to visit all nei...
Definition: intersectioniterator.hh:28
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:81
Intersection dereference() const
dereferencing
Definition: intersectioniterator.hh:105
void increment()
prefix increment
Definition: intersectioniterator.hh:93
bool equals(const MMeshInterfaceGridLeafIntersectionIterator &other) const
returns if iterators reference same intersection
Definition: intersectioniterator.hh:87
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:36
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 (Sep 4, 22:38, 2025)