5#ifndef DUNE_GEOGRID_INTERSECTION_HH
6#define DUNE_GEOGRID_INTERSECTION_HH
8#include <dune/grid/geometrygrid/declaration.hh>
9#include <dune/grid/geometrygrid/cornerstorage.hh>
20 template<
class Gr
id,
class HostIntersection >
23 typedef typename HostIntersection::Geometry HostGeometry;
24 typedef typename HostIntersection::LocalGeometry HostLocalGeometry;
26 typedef typename std::remove_const< Grid >::type::Traits Traits;
29 typedef typename Traits::ctype ctype;
31 static const int dimension = Traits::dimension;
32 static const int dimensionworld = Traits::dimensionworld;
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;
38 typedef typename Traits::template Codim< 0 >::Geometry ElementGeometry;
41 typedef GeoGrid::IntersectionCoordVector< Grid > CoordVector;
43 typedef typename Traits::template Codim< 0 >::EntityImpl EntityImpl;
45 typedef typename Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
46 typedef typename Traits::template Codim< 0 >::GeometryImpl ElementGeometryImpl;
53 explicit Intersection (
const HostIntersection &hostIntersection,
const ElementGeometryImpl &insideGeo )
54 : hostIntersection_( hostIntersection )
55 , insideGeo_ ( insideGeo )
59 explicit Intersection ( HostIntersection&& hostIntersection,
const ElementGeometryImpl &insideGeo )
60 : hostIntersection_(
std::move( hostIntersection ) )
61 , insideGeo_ ( insideGeo )
65 bool equals (
const Intersection &other)
const
67 return hostIntersection_ == other.hostIntersection_;
70 explicit operator bool ()
const {
return bool( hostIntersection_ ); }
72 Entity inside ()
const
74 return EntityImpl( insideGeo_, hostIntersection().inside() );
77 Entity outside ()
const
79 return EntityImpl( grid(), hostIntersection().outside() );
82 bool boundary ()
const {
return hostIntersection().boundary(); }
84 bool conforming ()
const {
return hostIntersection().conforming(); }
86 bool neighbor ()
const {
return hostIntersection().neighbor(); }
88 size_t boundarySegmentIndex ()
const
90 return hostIntersection().boundarySegmentIndex();
93 LocalGeometry geometryInInside ()
const
95 return hostIntersection().geometryInInside();
98 LocalGeometry geometryInOutside ()
const
100 return hostIntersection().geometryInOutside();
103 Geometry geometry ()
const
107 CoordVector coords( insideGeo_, geometryInInside() );
108 geo_ = GeometryImpl( grid(), type(), coords );
110 return Geometry( geo_ );
113 GeometryType type ()
const {
return hostIntersection().type(); }
115 int indexInInside ()
const
117 return hostIntersection().indexInInside();
120 int indexInOutside ()
const
122 return hostIntersection().indexInOutside();
125 FieldVector< ctype, dimensionworld >
126 integrationOuterNormal (
const FieldVector< ctype, dimension-1 > &local )
const
128 const LocalGeometry geoInInside = geometryInInside();
129 const int idxInInside = indexInInside();
131 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
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 );
137 FieldVector< ctype, dimensionworld > normal;
138 jit.mv( refNormal, normal );
140 normal *= geoInInside.volume() / refElement.template geometry< 1 >( idxInInside ).volume();
141 normal *= jit.detInv();
146 FieldVector< ctype, dimensionworld >
147 outerNormal (
const FieldVector< ctype, dimension-1 > &local )
const
149 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
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() );
155 FieldVector< ctype, dimensionworld > normal;
156 jit.mv( refNormal, normal );
160 FieldVector< ctype, dimensionworld >
161 unitOuterNormal (
const FieldVector< ctype, dimension-1 > &local )
const
163 FieldVector< ctype, dimensionworld > normal = outerNormal( local );
164 normal *= (ctype( 1 ) / normal.two_norm());
168 FieldVector< ctype, dimensionworld > centerUnitOuterNormal ()
const
171 return unitOuterNormal( refFace.position( 0, 0 ) );
174 const HostIntersection &hostIntersection ()
const
176 return hostIntersection_;
179 const Grid &grid ()
const {
return insideGeo_.grid(); }
182 HostIntersection hostIntersection_;
183 ElementGeometryImpl insideGeo_;
184 mutable GeometryImpl geo_;
@ 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:402
Dune namespace.
Definition: alignedallocator.hh:13