DUNE PDELab (2.8)

Dune::Intersection< GridImp, IntersectionImp > Class Template Reference

Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the domain boundary. More...

#include <dune/grid/common/intersection.hh>

Public Types

enum  { mydimension =GridImp::dimension-1 }
 Export dimension of the intersection. More...
 
enum  { dimensionworld =GridImp::dimensionworld }
 Export dimension of world. More...
 
typedef IntersectionImp Implementation
 type of underlying implementation More...
 
typedef GridImp::template Codim< 0 >::Entity Entity
 Type of entity that this Intersection belongs to.
 
typedef GridImp::template Codim< 1 >::Geometry Geometry
 Codim 1 geometry returned by geometry()
 
typedef Geometry::LocalCoordinate LocalCoordinate
 Type for vectors of coordinates on the intersection.
 
typedef Geometry::GlobalCoordinate GlobalCoordinate
 Type for normal vectors.
 
typedef GridImp::template Codim< 1 >::LocalGeometry LocalGeometry
 Codim 1 geometry returned by geometryInInside() and geometryInOutside()
 
typedef GridImp::ctype ctype
 Type of individual coefficients of coordinate vectors.
 

Public Member Functions

Implementationimpl ()
 access to the underlying implementation More...
 
const Implementationimpl () const
 access to the underlying implementation More...
 
bool boundary () const
 Return true if intersection is with interior or exterior boundary (see the cases above)
 
size_t boundarySegmentIndex () const
 index of the boundary segment within the macro grid More...
 
bool neighbor () const
 return true if intersection is shared with another element.
 
Entity inside () const
 return Entity on the inside of this intersection. That is the Entity where we started this.
 
Entity outside () const
 return Entity on the outside of this intersection. That is the neighboring Entity. More...
 
bool conforming () const
 Return true if intersection is conforming.
 
LocalGeometry geometryInInside () const
 geometrical information about this intersection in local coordinates of the inside() entity. More...
 
LocalGeometry geometryInOutside () const
 geometrical information about this intersection in local coordinates of the outside() entity. More...
 
Geometry geometry () const
 geometrical information about the intersection in global coordinates. More...
 
GeometryType type () const
 obtain the type of reference element for this intersection
 
int indexInInside () const
 Local index of codim 1 entity in the inside() entity where intersection is contained in. More...
 
int indexInOutside () const
 Local index of codim 1 entity in outside() entity where intersection is contained in. More...
 
GlobalCoordinate outerNormal (const LocalCoordinate &local) const
 Return an outer normal (length not necessarily 1) More...
 
GlobalCoordinate integrationOuterNormal (const LocalCoordinate &local) const
 return unit outer normal scaled with the integration element More...
 
GlobalCoordinate unitOuterNormal (const LocalCoordinate &local) const
 Return unit outer normal (length == 1) More...
 
GlobalCoordinate centerUnitOuterNormal () const
 Return unit outer normal (length == 1) More...
 
bool operator== (const Intersection &other) const
 Compares two intersections for equality.
 
bool operator!= (const Intersection &other) const
 Compares two intersections for inequality.
 
 Intersection ()
 Default constructor.
 
 Intersection (const Intersection &other)
 Copy constructor from an existing intersection.
 
 Intersection (Intersection &&other)
 Move constructor from an existing intersection.
 
Intersectionoperator= (const Intersection &other)
 Copy assignment operator from an existing intersection.
 
Intersectionoperator= (Intersection &&other)
 Move assignment operator from an existing intersection.
 

Implementor interface

 Intersection (const Implementation &impl)
 
 Intersection (Implementation &&impl)
 

Detailed Description

template<class GridImp, class IntersectionImp>
class Dune::Intersection< GridImp, IntersectionImp >

Intersection of a mesh entity of codimension 0 ("element") with a "neighboring" element or with the domain boundary.

Template Parameters
GridImpType that is a model of Dune::Grid
IntersectionImpClass template that is a model of Dune::Intersection

Overview

Intersections are codimension 1 objects. These intersections are accessed via an Intersection. This allows the implementation of non-matching grids, as one facet can now consist of several intersections. In a conforming mesh such an intersection corresponds to an entity of codimension 1 but in the general non-conforming case there will be no entity in the mesh that directly corresponds to the intersection. Thus, the Intersection describes these intersections implicitly.

Engine Concept

The Intersection class template wraps an object of type IntersectionImp and forwards all member function calls to corresponding members of this class. In that sense Intersection defines the interface and IntersectionImp supplies the implementation.

Methods neighbor and boundary

The following holds for both the level and the leaf intersection : The intersection is started on a codimension 0 entity of the grid. If this entity belongs to the interior or the overlap region (see. ???) then the union of all intersections is identical to the boundary of the entity. On ghost elements the only stops on the border of the domain, i.e., only on the intersections with entities in the interior region. Depending on the boolean values returned by the methods boundary() and neighbor() one can detect the position of an intersection relative to the boundary. In general the method boundary() returns true if and only if the intersection is part of the physical boundary of the domain. The method neighbor() returns true only if the method outside() has a well defined return value.

The following cases are possible if the intersection is started on an entity in the interior or overlap region. More details are given below:

intersectionneighbor()boundary()outside()
1with inner, overlap
or ghost entity
truefalse the neighbor entity
2on domain boundary falsetrueundefined
3on periodic boundary truetrueGhost-/Overlap cell
(with transformed geometry)
4on processor boundary false if grid has no ghosts
true otherwise
false ghost entity (if it exists)
Inner Intersections:
The type of the neighboring entity can be determined through methods defined on the outside entity.
Handling physical boundaries:
Data (like the type of boundary) can be attached to physical boundaries either using global coordinates or the intersection's boundary segment index.
The boundary segment indices form a local, zero-based, contiguous set of integer values. Each boundary segment on the macro level is assigned a unique index from this set, which is inherited by child boundary segments. The size of the boundary segment index set (i.e., the number of boundary indices on the macro level) can be determined through the method Grid::numBoundarySegments.
Note that the boundary segment index is not persistent over dynamic load balancing.
Handling periodic boundaries:
  • The Intersection stops at periodic boundaries
  • periodic grids are handled in correspondence to parallel grids
  • At the periodic boundary one can adjust an overlap- or ghost-layer.
  • outer() returns a ghost or overlap cell (for ghost and overlap look into the documentation of the parallel grid interface)
  • outer() cell has a periodically transformed geometry (so that one does not see a jump or something like that)
  • outer() cell has its own index
  • outer() cell has the same id as the corresponding "original" cell
Processor boundaries:
At processor boundaries, i.e. when an element has an intersection with another element in the sequential grid but this element is only stored in other processors the intersection stops but neither neighbor() nor boundary() are true.

Geometry of an intersection

The method geometry returns a geometry mapping the intersection as a codim one structure to global coordinates. The methods geometryInInside and geometryInOutside return geometries mapping the intersection into the reference elements of the originating entity and the neighboring entity, respectively. The indexInInside and indexInOutside methods return the codim one subentities which contain the intersection.

Member Typedef Documentation

◆ Implementation

template<class GridImp , class IntersectionImp >
typedef IntersectionImp Dune::Intersection< GridImp, IntersectionImp >::Implementation

type of underlying implementation

Warning
Implementation details may change without prior notification.

Member Enumeration Documentation

◆ anonymous enum

template<class GridImp , class IntersectionImp >
anonymous enum

Export dimension of the intersection.

Enumerator
mydimension 

intersection's dimension

◆ anonymous enum

template<class GridImp , class IntersectionImp >
anonymous enum

Export dimension of world.

Enumerator
dimensionworld 

dimension of world

Constructor & Destructor Documentation

◆ Intersection() [1/2]

template<class GridImp , class IntersectionImp >
Dune::Intersection< GridImp, IntersectionImp >::Intersection ( const Implementation impl)
inline

Copy Constructor from IntersectionImp

◆ Intersection() [2/2]

template<class GridImp , class IntersectionImp >
Dune::Intersection< GridImp, IntersectionImp >::Intersection ( Implementation &&  impl)
inline

Move Constructor from IntersectionImp

Member Function Documentation

◆ boundarySegmentIndex()

template<class GridImp , class IntersectionImp >
size_t Dune::Intersection< GridImp, IntersectionImp >::boundarySegmentIndex ( ) const
inline

index of the boundary segment within the macro grid

In many applications, special data needs to be attached to the boundary segments of the macro grid (e.g., a function selecting the boundary condition). Usually, this data is inherited by the children of the boundary segment.

In the DUNE framework, data is stored in arrays, addressed by an index, in this case the boundarySegmentIndex. The size of these arrays can be obtained by the Grid::numBoundarySegments method.

The indices returned by this method are consecutive, zero based, and local to the processor. Notice that these indices do not necessarily coincide with the insertion indices of the corresponding boundary segments as provided by the GridFactory.

Referenced by Dune::GridPtr< GridType >::parameters().

◆ centerUnitOuterNormal()

template<class GridImp , class IntersectionImp >
GlobalCoordinate Dune::Intersection< GridImp, IntersectionImp >::centerUnitOuterNormal ( ) const
inline

Return unit outer normal (length == 1)

The returned vector is the normal at the center() of the intersection's geometry. It is scaled to have unit length.

◆ geometry()

template<class GridImp , class IntersectionImp >
Geometry Dune::Intersection< GridImp, IntersectionImp >::geometry ( ) const
inline

geometrical information about the intersection in global coordinates.

This method returns a Geometry object that provides a mapping from local coordinates of the intersection to global (world) coordinates.

Note
If the returned geometry has type none then only a limited set of features is available for the geometry, i.e. center and volume.
Previously, the geometry was encapsulated in the intersection object and a const reference was returned.
The returned geometry object is guaranteed to remain valid until the grid is modified (or deleted).

◆ geometryInInside()

template<class GridImp , class IntersectionImp >
LocalGeometry Dune::Intersection< GridImp, IntersectionImp >::geometryInInside ( ) const
inline

geometrical information about this intersection in local coordinates of the inside() entity.

This method returns a Geometry object that provides a mapping from local coordinates of the intersection to local coordinates of the inside() entity.

Note
Previously, the geometry was encapsulated in the intersection object and a const reference was returned.
The returned geometry object is guaranteed to remain valid until the grid is modified (or deleted).

◆ geometryInOutside()

template<class GridImp , class IntersectionImp >
LocalGeometry Dune::Intersection< GridImp, IntersectionImp >::geometryInOutside ( ) const
inline

geometrical information about this intersection in local coordinates of the outside() entity.

This method returns a Geometry object that provides a mapping from local coordinates of the intersection to local coordinates of the outside() entity.

Note
Previously, the geometry was encapsulated in the intersection object and a const reference was returned.
The returned geometry object is guaranteed to remain valid until the grid is modified (or deleted).

◆ impl() [1/2]

template<class GridImp , class IntersectionImp >
Implementation & Dune::Intersection< GridImp, IntersectionImp >::impl ( )
inline

access to the underlying implementation

Warning
Implementation details may change without prior notification.

Referenced by Dune::YaspIntersectionIterator< GridImp >::dereference(), and Dune::YaspIntersectionIterator< GridImp >::increment().

◆ impl() [2/2]

template<class GridImp , class IntersectionImp >
const Implementation & Dune::Intersection< GridImp, IntersectionImp >::impl ( ) const
inline

access to the underlying implementation

Warning
Implementation details may change without prior notification.

◆ indexInInside()

template<class GridImp , class IntersectionImp >
int Dune::Intersection< GridImp, IntersectionImp >::indexInInside ( ) const
inline

Local index of codim 1 entity in the inside() entity where intersection is contained in.

This method returns the facet number with respect to the generic reference element.

Note
This index can be used with the inside entity's subEntity<1>() method to obtain the facet.
Returns
the index of the inside entity's facet containing this intersection (with respect to the generic reference element)

Referenced by Dune::Functions::SubEntityDOFs< GridView >::bind().

◆ indexInOutside()

template<class GridImp , class IntersectionImp >
int Dune::Intersection< GridImp, IntersectionImp >::indexInOutside ( ) const
inline

Local index of codim 1 entity in outside() entity where intersection is contained in.

This method returns the facet number with respect to the generic reference element.

Note
This index can be used with the outside entity's subEntity<1>() method to obtain the facet.
Returns
the index of the outside entity's facet containing this intersection (with respect to the generic reference element)

◆ integrationOuterNormal()

template<class GridImp , class IntersectionImp >
GlobalCoordinate Dune::Intersection< GridImp, IntersectionImp >::integrationOuterNormal ( const LocalCoordinate local) const
inline

return unit outer normal scaled with the integration element

The normal is scaled with the integration element of the intersection. This method is redundant but it may be more efficient to use this function rather than computing the integration element via geometry().

The returned vector may depend on local position within the intersection.

◆ outerNormal()

template<class GridImp , class IntersectionImp >
GlobalCoordinate Dune::Intersection< GridImp, IntersectionImp >::outerNormal ( const LocalCoordinate local) const
inline

Return an outer normal (length not necessarily 1)

The returned vector may depend on local position within the intersection.

◆ outside()

template<class GridImp , class IntersectionImp >
Entity Dune::Intersection< GridImp, IntersectionImp >::outside ( ) const
inline

return Entity on the outside of this intersection. That is the neighboring Entity.

Warning
Don't call this method if there is no neighboring Entity (neighbor() returns false). In this case the result is undefined.

◆ unitOuterNormal()

template<class GridImp , class IntersectionImp >
GlobalCoordinate Dune::Intersection< GridImp, IntersectionImp >::unitOuterNormal ( const LocalCoordinate local) const
inline

Return unit outer normal (length == 1)

The returned vector may depend on the local position within the intersection. It is scaled to have unit length.


The documentation for this class was generated from the following files:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 21, 23:30, 2024)