DUNE PDELab (git)

intersection.hh
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_GEOGRID_INTERSECTION_HH
6#define DUNE_GEOGRID_INTERSECTION_HH
7
8#include <dune/grid/geometrygrid/declaration.hh>
9#include <dune/grid/geometrygrid/cornerstorage.hh>
10
11namespace Dune
12{
13
14 namespace GeoGrid
15 {
16
17 // Intersection
18 // ------------
19
20 template< class Grid, class HostIntersection >
21 class Intersection
22 {
23 typedef typename HostIntersection::Geometry HostGeometry;
24 typedef typename HostIntersection::LocalGeometry HostLocalGeometry;
25
26 typedef typename std::remove_const< Grid >::type::Traits Traits;
27
28 public:
29 typedef typename Traits::ctype ctype;
30
31 static const int dimension = Traits::dimension;
32 static const int dimensionworld = Traits::dimensionworld;
33
34 typedef typename Traits::template Codim< 0 >::Entity Entity;
35 typedef typename Traits::template Codim< 1 >::Geometry Geometry;
36 typedef typename Traits::template Codim< 1 >::LocalGeometry LocalGeometry;
37
38 typedef typename Traits::template Codim< 0 >::Geometry ElementGeometry;
39
40 private:
41 typedef GeoGrid::IntersectionCoordVector< Grid > CoordVector;
42
43 typedef typename Traits::template Codim< 0 >::EntityImpl EntityImpl;
44
45 typedef typename Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
46 typedef typename Traits::template Codim< 0 >::GeometryImpl ElementGeometryImpl;
47
48 public:
49
50 Intersection()
51 {}
52
53 explicit Intersection ( const HostIntersection &hostIntersection, const ElementGeometryImpl &insideGeo )
54 : hostIntersection_( hostIntersection )
55 , insideGeo_ ( insideGeo )
56 , geo_( grid() )
57 {}
58
59 explicit Intersection ( HostIntersection&& hostIntersection, const ElementGeometryImpl &insideGeo )
60 : hostIntersection_( std::move( hostIntersection ) )
61 , insideGeo_ ( insideGeo )
62 , geo_( grid() )
63 {}
64
65 bool equals ( const Intersection &other) const
66 {
67 return hostIntersection_ == other.hostIntersection_;
68 }
69
70 explicit operator bool () const { return bool( hostIntersection_ ); }
71
72 Entity inside () const
73 {
74 return EntityImpl( insideGeo_, hostIntersection().inside() );
75 }
76
77 Entity outside () const
78 {
79 return EntityImpl( grid(), hostIntersection().outside() );
80 }
81
82 bool boundary () const { return hostIntersection().boundary(); }
83
84 bool conforming () const { return hostIntersection().conforming(); }
85
86 bool neighbor () const { return hostIntersection().neighbor(); }
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 LocalGeometry geoInInside = geometryInInside();
129 const int idxInInside = indexInInside();
130
131 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
132
133 FieldVector< ctype, dimension > x( geoInInside.global( local ) );
134 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
135 FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( idxInInside );
136
137 FieldVector< ctype, dimensionworld > normal;
138 jit.mv( refNormal, normal );
139 if( !conforming() )
140 normal *= geoInInside.volume() / refElement.template geometry< 1 >( idxInInside ).volume();
141 normal *= jit.detInv();
142 //normal *= insideGeo_.integrationElement( x );
143 return normal;
144 }
145
146 FieldVector< ctype, dimensionworld >
147 outerNormal ( const FieldVector< ctype, dimension-1 > &local ) const
148 {
149 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
150
151 FieldVector< ctype, dimension > x( geometryInInside().global( local ) );
152 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
153 FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( indexInInside() );
154
155 FieldVector< ctype, dimensionworld > normal;
156 jit.mv( refNormal, normal );
157 return normal;
158 }
159
160 FieldVector< ctype, dimensionworld >
161 unitOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
162 {
163 FieldVector< ctype, dimensionworld > normal = outerNormal( local );
164 normal *= (ctype( 1 ) / normal.two_norm());
165 return normal;
166 }
167
168 FieldVector< ctype, dimensionworld > centerUnitOuterNormal () const
169 {
170 auto refFace = referenceElement< ctype, dimension-1 >( type() );
171 return unitOuterNormal( refFace.position( 0, 0 ) );
172 }
173
174 const HostIntersection &hostIntersection () const
175 {
176 return hostIntersection_;
177 }
178
179 const Grid &grid () const { return insideGeo_.grid(); }
180
181 private:
182 HostIntersection hostIntersection_;
183 ElementGeometryImpl insideGeo_;
184 mutable GeometryImpl geo_;
185 };
186
187 } // namespace GeoGrid
188
189} // namespace Dune
190
191#endif // #ifndef DUNE_GEOGRID_INTERSECTION_HH
@ conforming
Output conforming data.
Definition: common.hh:73
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
unspecified value type referenceElement(T &&... t)
Returns a reference element for the objects t....
constexpr auto equals(T1 &&t1, T2 &&t2)
Equality comparison.
Definition: hybridutilities.hh:587
Dune namespace.
Definition: alignedallocator.hh:13
STL namespace.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)