3#ifndef DUNE_GEOGRID_INTERSECTION_HH
4#define DUNE_GEOGRID_INTERSECTION_HH
6#include <dune/grid/geometrygrid/declaration.hh>
7#include <dune/grid/geometrygrid/cornerstorage.hh>
18 template<
class Gr
id,
class HostIntersection >
21 typedef typename HostIntersection::Geometry HostGeometry;
22 typedef typename HostIntersection::LocalGeometry HostLocalGeometry;
24 typedef typename std::remove_const< Grid >::type::Traits Traits;
27 typedef typename Traits::ctype ctype;
29 static const int dimension = Traits::dimension;
30 static const int dimensionworld = Traits::dimensionworld;
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;
36 typedef typename Traits::template Codim< 0 >::Geometry ElementGeometry;
39 typedef GeoGrid::IntersectionCoordVector< Grid > CoordVector;
41 typedef typename Traits::template Codim< 0 >::EntityImpl EntityImpl;
43 typedef typename Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
44 typedef typename Traits::template Codim< 0 >::GeometryImpl ElementGeometryImpl;
51 explicit Intersection (
const HostIntersection &hostIntersection,
const ElementGeometryImpl &insideGeo )
52 : hostIntersection_( hostIntersection )
53 , insideGeo_ ( insideGeo )
57 explicit Intersection ( HostIntersection&& hostIntersection,
const ElementGeometryImpl &insideGeo )
58 : hostIntersection_(
std::move( hostIntersection ) )
59 , insideGeo_ ( insideGeo )
63 bool equals (
const Intersection &other)
const
65 return hostIntersection_ == other.hostIntersection_;
68 explicit operator bool ()
const {
return bool( hostIntersection_ ); }
70 Entity inside ()
const
72 return EntityImpl( insideGeo_, hostIntersection().inside() );
75 Entity outside ()
const
77 return EntityImpl( grid(), hostIntersection().outside() );
80 bool boundary ()
const {
return hostIntersection().boundary(); }
82 bool conforming ()
const {
return hostIntersection().conforming(); }
84 bool neighbor ()
const {
return hostIntersection().neighbor(); }
86 size_t boundarySegmentIndex ()
const
88 return hostIntersection().boundarySegmentIndex();
91 LocalGeometry geometryInInside ()
const
93 return hostIntersection().geometryInInside();
96 LocalGeometry geometryInOutside ()
const
98 return hostIntersection().geometryInOutside();
101 Geometry geometry ()
const
105 CoordVector coords( insideGeo_, geometryInInside() );
106 geo_ = GeometryImpl( grid(), type(), coords );
108 return Geometry( geo_ );
111 GeometryType type ()
const {
return hostIntersection().type(); }
113 int indexInInside ()
const
115 return hostIntersection().indexInInside();
118 int indexInOutside ()
const
120 return hostIntersection().indexInOutside();
123 FieldVector< ctype, dimensionworld >
124 integrationOuterNormal (
const FieldVector< ctype, dimension-1 > &local )
const
126 const LocalGeometry geoInInside = geometryInInside();
127 const int idxInInside = indexInInside();
129 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
131 FieldVector< ctype, dimension > x( geoInInside.global( local ) );
132 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
133 FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( idxInInside );
135 FieldVector< ctype, dimensionworld > normal;
136 jit.mv( refNormal, normal );
138 normal *= geoInInside.volume() / refElement.template geometry< 1 >( idxInInside ).volume();
139 normal *= jit.detInv();
144 FieldVector< ctype, dimensionworld >
145 outerNormal (
const FieldVector< ctype, dimension-1 > &local )
const
147 auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
149 FieldVector< ctype, dimension > x( geometryInInside().global( local ) );
150 const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
151 FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( indexInInside() );
153 FieldVector< ctype, dimensionworld > normal;
154 jit.mv( refNormal, normal );
158 FieldVector< ctype, dimensionworld >
159 unitOuterNormal (
const FieldVector< ctype, dimension-1 > &local )
const
161 FieldVector< ctype, dimensionworld > normal = outerNormal( local );
162 normal *= (ctype( 1 ) / normal.two_norm());
166 FieldVector< ctype, dimensionworld > centerUnitOuterNormal ()
const
169 return unitOuterNormal( refFace.position( 0, 0 ) );
172 const HostIntersection &hostIntersection ()
const
174 return hostIntersection_;
177 const Grid &grid ()
const {
return insideGeo_.grid(); }
180 HostIntersection hostIntersection_;
181 ElementGeometryImpl insideGeo_;
182 mutable GeometryImpl geo_;
@ conforming
Output conforming data.
Definition: common.hh:72
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:180
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:401
Dune namespace.
Definition: alignedallocator.hh:14