DUNE-FEM (unstable)

deaditerator.hh
1#ifndef DUNE_FEM_GRIDPART_COMMON_DEADINTERSECTIONITERATOR_HH
2#define DUNE_FEM_GRIDPART_COMMON_DEADINTERSECTIONITERATOR_HH
3
4#include <type_traits>
5
6#include <dune/grid/common/entityiterator.hh>
7#include <dune/grid/common/intersection.hh>
8#include <dune/grid/common/intersectioniterator.hh>
9
10namespace Dune
11{
12
13 namespace Fem
14 {
15
16 // DeadIterator
17 // ------------
18
19 template< class E >
20 struct DeadIterator
21 {
22 typedef E Entity;
23
24 static const int codimension = Entity::codimension;
25
26 void increment ()
27 {
28 DUNE_THROW( InvalidStateException, "Trying to increment a dead iterator." );
29 }
30
31 Entity dereference () const
32 {
33 DUNE_THROW( InvalidStateException, "Trying to dereference a dead iterator." );
34 }
35
36 bool equals ( const DeadIterator & ) const
37 {
38 DUNE_THROW( InvalidStateException, "Trying to compare a dead iterator." );
39 }
40 };
41
42
43
44 // DeadIntersection
45 // ----------------
46
47 template< class GridFamily >
48 class DeadIntersection
49 {
50 typedef typename std::remove_const< GridFamily >::type::Traits Traits;
51
52 public:
53 typedef typename std::remove_const< GridFamily >::type::ctype ctype;
54
55 static const int dimension = std::remove_const< GridFamily >::type::dimension;
56 static const int dimensionworld = std::remove_const< GridFamily >::type::dimensionworld;
57
58 typedef typename Traits::template Codim< 0 >::Entity Entity;
59 typedef typename Traits::template Codim< 1 >::Geometry Geometry;
60 typedef typename Traits::template Codim< 1 >::LocalGeometry LocalGeometry;
61
62 Entity inside () const
63 {
64 DUNE_THROW( InvalidStateException, "Call to inside on dead intersection." );
65 }
66
67 Entity outside () const
68 {
69 DUNE_THROW( InvalidStateException, "Call to outside on dead intersection." );
70 }
71
72 bool boundary () const
73 {
74 DUNE_THROW( InvalidStateException, "Call to boundary on dead intersection." );
75 }
76
77 bool conforming () const
78 {
79 DUNE_THROW( InvalidStateException, "Call to conforming on dead intersection." );
80 }
81
82 bool neighbor () const
83 {
84 DUNE_THROW( InvalidStateException, "Call to neighbor on dead intersection." );
85 }
86
87 int boundaryId () const
88 {
89 DUNE_THROW( InvalidStateException, "Call to boundaryId on dead intersection." );
90 }
91
92 size_t boundarySegmentIndex () const
93 {
94 DUNE_THROW( InvalidStateException, "Call to boundarySegmentIndex on dead intersection." );
95 }
96
97 const LocalGeometry &geometryInInside () const
98 {
99 DUNE_THROW( InvalidStateException, "Call to geometryInInside on dead intersection." );
100 }
101
102 const LocalGeometry &geometryInOutside () const
103 {
104 DUNE_THROW( InvalidStateException, "Call to geometryInOutside on dead intersection." );
105 }
106
107 const Geometry &geometry () const
108 {
109 DUNE_THROW( InvalidStateException, "Call to geometry on dead intersection." );
110 }
111
112 GeometryType type () const
113 {
114 DUNE_THROW( InvalidStateException, "Call to type on dead intersection." );
115 }
116
117 int indexInInside () const
118 {
119 DUNE_THROW( InvalidStateException, "Call to indexInInside on dead intersection." );
120 }
121
122 int indexInOutside () const
123 {
124 DUNE_THROW( InvalidStateException, "Call to indexInOutside on dead intersection." );
125 }
126
127 FieldVector< ctype, dimensionworld >
128 integrationOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
129 {
130 DUNE_THROW( InvalidStateException, "Call to integrationOuterNormal on dead intersection." );
131 }
132
133 FieldVector< ctype, dimensionworld >
134 outerNormal ( const FieldVector< ctype, dimension-1 > &local ) const
135 {
136 DUNE_THROW( InvalidStateException, "Call to outerNormal on dead intersection." );
137 }
138
139 FieldVector< ctype, dimensionworld >
140 unitOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
141 {
142 DUNE_THROW( InvalidStateException, "Call to unitOuterNormal on dead intersection." );
143 }
144
145 FieldVector< ctype, dimensionworld > centerUnitOuterNormal () const
146 {
147 DUNE_THROW( InvalidStateException, "Call to centerUnitOuterNormal on dead intersection." );
148 }
149 };
150
151
152
153 // DeadIntersectionIterator
154 // ------------------------
155
156 template< class GridFamily >
157 class DeadIntersectionIterator
158 {
159 typedef DeadIntersectionIterator< GridFamily > ThisType;
160
161 typedef DeadIntersection< const GridFamily > DeadIntersectionType;
162
163 public:
165
166 DeadIntersectionIterator ()
167 {}
168
169 bool equals ( const ThisType &other ) const
170 {
171 return true;
172 }
173
174 void increment ()
175 {
176 DUNE_THROW( InvalidStateException, "Trying to increment a dead intersection iterator." );
177 }
178
179 const Intersection &dereference () const
180 {
181 DUNE_THROW( InvalidStateException, "Trying to dereference a dead intersection iterator." );
182 }
183 };
184
185 } // namespace Fem
186
187} // namespace Dune
188
189#endif // #ifndef DUNE_FEM_GRIDPART_COMMON_DEADINTERSECTIONITERATOR_HH
static constexpr int codimension
Know your own codimension.
Definition: entity.hh:106
Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the d...
Definition: intersection.hh:164
@ conforming
Output conforming data.
Definition: common.hh:73
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto equals
Function object for performing equality comparison.
Definition: hybridutilities.hh:572
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 27, 22:29, 2024)