1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191
|
// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifndef DUNE_GEOGRID_INTERSECTION_HH
#define DUNE_GEOGRID_INTERSECTION_HH
#include <dune/grid/geometrygrid/declaration.hh>
#include <dune/grid/geometrygrid/cornerstorage.hh>
namespace Dune
{
namespace GeoGrid
{
// Intersection
// ------------
template< class Grid, class HostIntersection >
class Intersection
{
typedef typename HostIntersection::Geometry HostGeometry;
typedef typename HostIntersection::LocalGeometry HostLocalGeometry;
typedef typename std::remove_const< Grid >::type::Traits Traits;
public:
typedef typename Traits::ctype ctype;
static const int dimension = Traits::dimension;
static const int dimensionworld = Traits::dimensionworld;
typedef typename Traits::template Codim< 0 >::Entity Entity;
typedef typename Traits::template Codim< 1 >::Geometry Geometry;
typedef typename Traits::template Codim< 1 >::LocalGeometry LocalGeometry;
typedef typename Traits::template Codim< 0 >::Geometry ElementGeometry;
private:
typedef GeoGrid::IntersectionCoordVector< Grid > CoordVector;
typedef typename Traits::template Codim< 0 >::EntityImpl EntityImpl;
typedef typename Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
typedef typename Traits::template Codim< 0 >::GeometryImpl ElementGeometryImpl;
public:
Intersection()
{}
explicit Intersection ( const HostIntersection &hostIntersection, const ElementGeometryImpl &insideGeo )
: hostIntersection_( hostIntersection )
, insideGeo_ ( insideGeo )
, geo_( grid() )
{}
explicit Intersection ( HostIntersection&& hostIntersection, const ElementGeometryImpl &insideGeo )
: hostIntersection_( std::move( hostIntersection ) )
, insideGeo_ ( insideGeo )
, geo_( grid() )
{}
bool equals ( const Intersection &other) const
{
return hostIntersection_ == other.hostIntersection_;
}
explicit operator bool () const { return bool( hostIntersection_ ); }
Entity inside () const
{
return EntityImpl( insideGeo_, hostIntersection().inside() );
}
Entity outside () const
{
return EntityImpl( grid(), hostIntersection().outside() );
}
bool boundary () const { return hostIntersection().boundary(); }
bool conforming () const { return hostIntersection().conforming(); }
bool neighbor () const { return hostIntersection().neighbor(); }
size_t boundarySegmentIndex () const
{
return hostIntersection().boundarySegmentIndex();
}
LocalGeometry geometryInInside () const
{
return hostIntersection().geometryInInside();
}
LocalGeometry geometryInOutside () const
{
return hostIntersection().geometryInOutside();
}
Geometry geometry () const
{
if( !geo_ )
{
CoordVector coords( insideGeo_, geometryInInside() );
geo_ = GeometryImpl( grid(), type(), coords );
}
return Geometry( geo_ );
}
GeometryType type () const { return hostIntersection().type(); }
int indexInInside () const
{
return hostIntersection().indexInInside();
}
int indexInOutside () const
{
return hostIntersection().indexInOutside();
}
FieldVector< ctype, dimensionworld >
integrationOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
{
const LocalGeometry geoInInside = geometryInInside();
const int idxInInside = indexInInside();
auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
FieldVector< ctype, dimension > x( geoInInside.global( local ) );
const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( idxInInside );
FieldVector< ctype, dimensionworld > normal;
jit.mv( refNormal, normal );
if( !conforming() )
normal *= geoInInside.volume() / refElement.template geometry< 1 >( idxInInside ).volume();
normal *= jit.detInv();
//normal *= insideGeo_.integrationElement( x );
return normal;
}
FieldVector< ctype, dimensionworld >
outerNormal ( const FieldVector< ctype, dimension-1 > &local ) const
{
auto refElement = referenceElement< ctype, dimension >( insideGeo_.type() );
FieldVector< ctype, dimension > x( geometryInInside().global( local ) );
const typename ElementGeometryImpl::JacobianInverseTransposed &jit = insideGeo_.jacobianInverseTransposed( x );
FieldVector< ctype, dimension > refNormal = refElement.integrationOuterNormal( indexInInside() );
FieldVector< ctype, dimensionworld > normal;
jit.mv( refNormal, normal );
return normal;
}
FieldVector< ctype, dimensionworld >
unitOuterNormal ( const FieldVector< ctype, dimension-1 > &local ) const
{
FieldVector< ctype, dimensionworld > normal = outerNormal( local );
normal *= (ctype( 1 ) / normal.two_norm());
return normal;
}
FieldVector< ctype, dimensionworld > centerUnitOuterNormal () const
{
auto refFace = referenceElement< ctype, dimension-1 >( type() );
return unitOuterNormal( refFace.position( 0, 0 ) );
}
const HostIntersection &hostIntersection () const
{
return hostIntersection_;
}
const Grid &grid () const { return insideGeo_.grid(); }
private:
HostIntersection hostIntersection_;
ElementGeometryImpl insideGeo_;
mutable GeometryImpl geo_;
};
} // namespace GeoGrid
} // namespace Dune
#endif // #ifndef DUNE_GEOGRID_INTERSECTION_HH
|