Dune Core Modules (2.4.1)

intersection.hh
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_GEOGRID_INTERSECTION_HH
4#define DUNE_GEOGRID_INTERSECTION_HH
5
6#include <dune/grid/geometrygrid/declaration.hh>
7#include <dune/grid/geometrygrid/cornerstorage.hh>
8
9namespace Dune
10{
11
12 namespace GeoGrid
13 {
14
15 // Intersection
16 // ------------
17
18 template< class Grid, class HostIntersection >
19 class Intersection
20 {
21 typedef typename HostIntersection::Geometry HostGeometry;
22 typedef typename HostIntersection::LocalGeometry HostLocalGeometry;
23
24 typedef typename remove_const< Grid >::type::Traits Traits;
25
26 public:
27 typedef typename Traits::ctype ctype;
28
29 static const int dimension = Traits::dimension;
30 static const int dimensionworld = Traits::dimensionworld;
31
32 typedef typename Traits::template Codim< 0 >::Entity Entity;
33 typedef typename Traits::template Codim< 1 >::Geometry Geometry;
34 typedef typename Traits::template Codim< 1 >::LocalGeometry LocalGeometry;
35
36 typedef typename Traits::template Codim< 0 >::Geometry ElementGeometry;
37
38 private:
39 typedef GeoGrid::IntersectionCoordVector< Grid > CoordVector;
40
41 typedef typename Traits::template Codim< 0 >::EntityImpl EntityImpl;
42
43 typedef typename Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
44 typedef typename Traits::template Codim< 0 >::GeometryImpl ElementGeometryImpl;
45
46 public:
47
48 Intersection()
49 {}
50
51 explicit Intersection ( const HostIntersection &hostIntersection, const ElementGeometryImpl &insideGeo )
52 : hostIntersection_( hostIntersection )
53 , insideGeo_ ( insideGeo )
54 , geo_( grid() )
55 {}
56
57 explicit Intersection ( HostIntersection&& hostIntersection, const ElementGeometryImpl &insideGeo )
58 : hostIntersection_( std::move( hostIntersection ) )
59 , insideGeo_ ( insideGeo )
60 , geo_( grid() )
61 {}
62
63 bool equals ( const Intersection &other) const
64 {
65 return hostIntersection_ == other.hostIntersection_;
66 }
67
68 operator bool () const { return bool( hostIntersection_ ); }
69
70 Entity inside () const
71 {
72 return EntityImpl( insideGeo_, hostIntersection().inside() );
73 }
74
75 Entity outside () const
76 {
77 return EntityImpl( grid(), hostIntersection().outside() );
78 }
79
80 bool boundary () const { return hostIntersection().boundary(); }
81
82 bool conforming () const { return hostIntersection().conforming(); }
83
84 bool neighbor () const { return hostIntersection().neighbor(); }
85
86 int boundaryId () const { return hostIntersection().boundaryId(); }
87
88 size_t boundarySegmentIndex () const
89 {
90 return hostIntersection().boundarySegmentIndex();
91 }
92
93 LocalGeometry geometryInInside () const
94 {
95 return hostIntersection().geometryInInside();
96 }
97
98 LocalGeometry geometryInOutside () const
99 {
100 return hostIntersection().geometryInOutside();
101 }
102
103 Geometry geometry () const
104 {
105 if( !geo_ )
106 {
107 CoordVector coords( insideGeo_, geometryInInside() );
108 geo_ = GeometryImpl( grid(), type(), coords );
109 }
110 return Geometry( geo_ );
111 }
112
113 GeometryType type () const { return hostIntersection().type(); }
114
115 int indexInInside () const
116 {
117 return hostIntersection().indexInInside();
118 }
119
120 int indexInOutside () const
121 {
122 return hostIntersection().indexInOutside();
123 }
124
125 FieldVector< ctype, dimensionworld >
126 integrationOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
127 {
128 const ReferenceElement< ctype, dimension > &refElement
130
131 FieldVector< ctype, dimension > x( geometryInInside().global( local ) );
132 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
133 const FieldVector< ctype, dimension > &refNormal = refElement.integrationOuterNormal( indexInInside() );
134
135 FieldVector< ctype, dimensionworld > normal;
136 jit.mv( refNormal, normal );
137 normal *= jit.detInv();
138 //normal *= insideGeo_.integrationElement( x );
139 return normal;
140 }
141
142 FieldVector< ctype, dimensionworld >
143 outerNormal ( const FieldVector< ctype, dimension-1 > &local ) const
144 {
145 const ReferenceElement< ctype, dimension > &refElement
147
148 FieldVector< ctype, dimension > x( geometryInInside().global( local ) );
149 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
150 const FieldVector< ctype, dimension > &refNormal = refElement.integrationOuterNormal( indexInInside() );
151
152 FieldVector< ctype, dimensionworld > normal;
153 jit.mv( refNormal, normal );
154 return normal;
155 }
156
157 FieldVector< ctype, dimensionworld >
158 unitOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
159 {
160 FieldVector< ctype, dimensionworld > normal = outerNormal( local );
161 normal *= (ctype( 1 ) / normal.two_norm());
162 return normal;
163 }
164
165 FieldVector< ctype, dimensionworld > centerUnitOuterNormal () const
166 {
167 const ReferenceElement< ctype, dimension-1 > &refFace
169 return unitOuterNormal( refFace.position( 0, 0 ) );
170 }
171
172 const HostIntersection &hostIntersection () const
173 {
174 return hostIntersection_;
175 }
176
177 const Grid &grid () const { return insideGeo_.grid(); }
178
179 private:
180 HostIntersection hostIntersection_;
181 ElementGeometryImpl insideGeo_;
182 mutable GeometryImpl geo_;
183 };
184
185 } // namespace GeoGrid
186
187} // namespace Dune
188
189#endif // #ifndef DUNE_GEOGRID_INTERSECTION_HH
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
Dune namespace.
Definition: alignment.hh:10
STL namespace.
static const ReferenceElement< ctype, dim > & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:484
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 22, 23:33, 2024)