dune-grid-glue  2.6-git
Namespaces | Classes | Typedefs | Functions
Dune::GridGlue Namespace Reference

Namespaces

 AreaWriterImplementation
 
 ProjectionImplementation
 
 ProjectionWriterImplementation
 

Classes

class  Codim0Extractor
 
class  Codim1Extractor
 
class  CommDataHandle
 describes the features of a data handle for communication in parallel runs using the GridGlue::communicate methods. More...
 
struct  CommInfo
 collects all GridGlue data requried for communication More...
 
class  CommunicationOperator
 forward gather scatter to user defined CommInfo class More...
 
class  ComputationMethod
 
class  ConformingMerge
 Implementation of the Merger concept for conforming interfaces. More...
 
class  ContactMerge
 Merge two codimension-1 surfaces that may be a positive distance apart. More...
 
class  Extractor
 Provides codimension-independent methods for grid extraction. More...
 
struct  GlobalId
 
class  GridGlue
 sequential adapter to couple two grids at specified close together boundaries More...
 
class  GridGlueAmiraWriter
 Write remote intersections to a AmiraMesh file for debugging purposes. More...
 
class  GridGlueVtkWriter
 Write remote intersections to a vtk file for debugging purposes. More...
 
class  Intersection
 The intersection of two entities of the two patches of a GridGlue. More...
 
class  IntersectionComputation
 Intersection computation method for two elements of arbitrary dimension. More...
 
class  IntersectionData
 storage class for Dune::GridGlue::Intersection related data More...
 
class  IntersectionIndexSet
 
class  IntersectionIterator
 
class  IntersectionList
 
class  IntersectionListProvider
 
struct  IntersectionTraits
 
class  Merger
 Abstract base for all classes that take extracted grids and build sets of intersections. More...
 
class  OverlappingMerge
 Computing overlapping grid intersections for grids of different dimensions. More...
 
class  Projection
 Projection of a line (triangle) on another line (triangle). More...
 
struct  Reverse
 
class  SimplexMethod
 
class  SimplexMethod< dimWorld, 0, 0, T >
 
class  SimplexMethod< dimWorld, 0, 1, T >
 
class  SimplexMethod< dimWorld, 0, 2, T >
 
class  SimplexMethod< dimWorld, 0, 3, T >
 
class  SimplexMethod< dimWorld, 1, 1, T >
 
class  SimplexMethod< dimWorld, 1, 2, T >
 
class  SimplexMethod< dimWorld, 1, 3, T >
 
class  SimplexMethod< dimWorld, 2, 2, T >
 
class  SimplexMethod< dimWorld, 2, 3, T >
 
class  SimplexMethod< dimWorld, 3, 3, T >
 
class  SimplicialIntersectionListProvider
 
class  StandardMerge
 Common base class for many merger implementations: produce pairs of entities that may intersect. More...
 
class  StreamingMessageBuffer
 
class  VtkSurfaceWriter
 

Typedefs

typedef std::pair< int, int > RankPair
 
typedef CommunicationOperator< Dune::ForwardCommunication > ForwardOperator
 
typedef CommunicationOperator< Dune::BackwardCommunication > BackwardOperator
 

Functions

template<typename T >
void printVector (const std::vector< T > &v, std::string name, int rank)
 
std::ostream & operator<< (std::ostream &os, const GlobalId &id)
 
template<... >
IteratorRange<... > intersections (const GridGlue<... > &glue, const Reverse<... > &reverse=!reversed)
 Iterate over all intersections of a GridGlue. More...
 
template<int side, typename Glue >
void write_glue_area_vtk (const Glue &glue, std::ostream &out)
 
template<int side, typename Glue >
void write_glue_area_vtk (const Glue &glue, const std::string &filename)
 
template<typename Glue >
void write_glue_areas_vtk (const Glue &glue, const std::string &base)
 
template<class T , int dim>
static Dune::FieldVector< T, dim > crossProduct (const Dune::FieldVector< T, dim > &a, const Dune::FieldVector< T, dim > &b)
 compute cross product More...
 
template<typename Coordinate , typename Corners , typename Normals >
void write (const Projection< Coordinate > &projection, const Corners &corners, const Normals &normals, std::ostream &out)
 write projection in VTK format More...
 
template<typename Coordinate , typename Corners , typename Normals >
void write (const Projection< Coordinate > &projection, const Corners &corners, const Normals &normals, const std::string &filename)
 write projection in VTK format More...
 
template<typename Coordinate , typename Corners , typename Normals >
void print (const Projection< Coordinate > &projection, const Corners &corners, const Normals &normals)
 Print information about the projection to std::cout stream. More...
 
template<class V >
int insertPoint (const V p, std::vector< V > &P)
 
template<int dimworld, typename T >
void simplexSubdivision (std::integral_constant< int, 0 >, const std::vector< Dune::FieldVector< T, dimworld > > &elementCorners, std::vector< std::vector< unsigned int > > &subElements, std::vector< std::vector< int > > &faceIds)
 
template<int dimworld, typename T >
void simplexSubdivision (std::integral_constant< int, 1 >, const std::vector< Dune::FieldVector< T, dimworld > > &elementCorners, std::vector< std::vector< unsigned int > > &subElements, std::vector< std::vector< int > > &faceIds)
 
template<int dimworld, typename T >
void simplexSubdivision (std::integral_constant< int, 2 >, const std::vector< Dune::FieldVector< T, dimworld > > &elementCorners, std::vector< std::vector< unsigned int > > &subElements, std::vector< std::vector< int > > &faceIds)
 
template<int dimworld, typename T >
void simplexSubdivision (std::integral_constant< int, 3 >, const std::vector< Dune::FieldVector< T, dimworld > > &elementCorners, std::vector< std::vector< unsigned int > > &subElements, std::vector< std::vector< int > > &faceIds)
 
 STANDARD_MERGE_INSTANTIATE (double, 1, 1, 1)
 
 STANDARD_MERGE_INSTANTIATE (double, 2, 2, 2)
 
 STANDARD_MERGE_INSTANTIATE (double, 3, 3, 3)
 

Typedef Documentation

◆ BackwardOperator

typedef CommunicationOperator<Dune::BackwardCommunication> Dune::GridGlue::BackwardOperator

◆ ForwardOperator

typedef CommunicationOperator<Dune::ForwardCommunication> Dune::GridGlue::ForwardOperator

◆ RankPair

typedef std::pair<int, int> Dune::GridGlue::RankPair

Function Documentation

◆ crossProduct()

template<class T , int dim>
static Dune::FieldVector<T,dim> Dune::GridGlue::crossProduct ( const Dune::FieldVector< T, dim > &  a,
const Dune::FieldVector< T, dim > &  b 
)
static

compute cross product

Returns
a × b

◆ insertPoint()

template<class V >
int Dune::GridGlue::insertPoint ( const V  p,
std::vector< V > &  P 
)
inline

◆ intersections()

template<... >
IteratorRange<... > intersections ( const GridGlue<... > &  glue,
const Reverse<... > &  reverse = !reversed 
)

Iterate over all intersections of a GridGlue.

This function returns an object representing the range of intersections with respect to the GridGlue glue. Its main purpose is to enable iteration over these intersections by means of a range-based for loop:

// Iterate over all intersections of a GridGlue in various ways
using Dune::GridGlue::reversed;
GridGlue<...> glue;
for (const auto& in : intersections(glue)) { ... }
for (const auto& in : intersections(glue, reversed)) { ... }
for (const auto& in : intersections(glue, !reversed)) { ... }
for (const auto& in : intersections(glue, Reversed<true>())) { ... }

The in- and outside of the intersection can be reversed by passing reversed as the second argument. The fourth form can be used in case a template parameter for reversal is required.

Since
dune-common 2.4
Parameters
glueGridGlue to obtain the intersections from
reverseTag to indicate reversal of in- and outside of intersections
Returns
an unspecified object that is guaranteed to fulfill the interface of IteratorRange and that can be iterated over using a range-based for loop.
See also
Dune::GridGlue::Intersection

◆ operator<<()

std::ostream& Dune::GridGlue::operator<< ( std::ostream &  os,
const GlobalId id 
)
inline

◆ print()

template<typename Coordinate , typename Corners , typename Normals >
void Dune::GridGlue::print ( const Projection< Coordinate > &  projection,
const Corners &  corners,
const Normals &  normals 
)

Print information about the projection to std::cout stream.

This method is mainly intended to be used for debugging.

Parameters
projectionprojection result
cornerscorners of the projected triangles
normalsnormals of the projected triangles

◆ printVector()

template<typename T >
void Dune::GridGlue::printVector ( const std::vector< T > &  v,
std::string  name,
int  rank 
)

◆ simplexSubdivision() [1/4]

template<int dimworld, typename T >
void Dune::GridGlue::simplexSubdivision ( std::integral_constant< int, 0 >  ,
const std::vector< Dune::FieldVector< T, dimworld > > &  elementCorners,
std::vector< std::vector< unsigned int > > &  subElements,
std::vector< std::vector< int > > &  faceIds 
)
inline

◆ simplexSubdivision() [2/4]

template<int dimworld, typename T >
void Dune::GridGlue::simplexSubdivision ( std::integral_constant< int, 1 >  ,
const std::vector< Dune::FieldVector< T, dimworld > > &  elementCorners,
std::vector< std::vector< unsigned int > > &  subElements,
std::vector< std::vector< int > > &  faceIds 
)
inline

◆ simplexSubdivision() [3/4]

template<int dimworld, typename T >
void Dune::GridGlue::simplexSubdivision ( std::integral_constant< int, 2 >  ,
const std::vector< Dune::FieldVector< T, dimworld > > &  elementCorners,
std::vector< std::vector< unsigned int > > &  subElements,
std::vector< std::vector< int > > &  faceIds 
)
inline

◆ simplexSubdivision() [4/4]

template<int dimworld, typename T >
void Dune::GridGlue::simplexSubdivision ( std::integral_constant< int, 3 >  ,
const std::vector< Dune::FieldVector< T, dimworld > > &  elementCorners,
std::vector< std::vector< unsigned int > > &  subElements,
std::vector< std::vector< int > > &  faceIds 
)
inline

◆ STANDARD_MERGE_INSTANTIATE() [1/3]

Dune::GridGlue::STANDARD_MERGE_INSTANTIATE ( double  ,
,
,
 
)

◆ STANDARD_MERGE_INSTANTIATE() [2/3]

Dune::GridGlue::STANDARD_MERGE_INSTANTIATE ( double  ,
,
,
 
)

◆ STANDARD_MERGE_INSTANTIATE() [3/3]

Dune::GridGlue::STANDARD_MERGE_INSTANTIATE ( double  ,
,
,
 
)

◆ write() [1/2]

template<typename Coordinate , typename Corners , typename Normals >
void Dune::GridGlue::write ( const Projection< Coordinate > &  projection,
const Corners &  corners,
const Normals &  normals,
const std::string &  filename 
)

write projection in VTK format

This function works the same as write(const Projection<Coordinate>&, const Corners&, const Normals&, std::ostream&) but creates a new file with the given name.

◆ write() [2/2]

template<typename Coordinate , typename Corners , typename Normals >
void Dune::GridGlue::write ( const Projection< Coordinate > &  projection,
const Corners &  corners,
const Normals &  normals,
std::ostream &  out 
)

write projection in VTK format

This file writes the projection information in VTK format to the given output stream. It is intended to be used for debugging so the details of the output might change at any time.

Note
This currently only works in 3D!
Parameters
projectionprojection result
cornerscorners of the projected triangles
normalsnormals of the projected triangles
outoutput stream

◆ write_glue_area_vtk() [1/2]

template<int side, typename Glue >
void Dune::GridGlue::write_glue_area_vtk ( const Glue &  glue,
const std::string &  filename 
)

◆ write_glue_area_vtk() [2/2]

template<int side, typename Glue >
void Dune::GridGlue::write_glue_area_vtk ( const Glue &  glue,
std::ostream &  out 
)

◆ write_glue_areas_vtk()

template<typename Glue >
void Dune::GridGlue::write_glue_areas_vtk ( const Glue &  glue,
const std::string &  base 
)

Write area covered by grid glue.

Write a VTK file for each side that indicate where the grid glue intersections are defined. A cell data field is provided that gives the relative area that is covered by glue intersections: if the value is zero no intersections are defined on the element, if the value is one the entire element should be covered (provided the glue is injective). Note that also values greater than one can be reached if the mapping is not injective.

This method is intended to be used for debugging purposes only.

Parameters
glueGridGlue for which the coverage by intersections should be computed
baseprefix of filenames to use. The generated files are named base-inside.vtk and base-outside.vtk.
Dune::GridGlue::intersections
IteratorRange<... > intersections(const GridGlue<... > &glue, const Reverse<... > &reverse=!reversed)
Iterate over all intersections of a GridGlue.
Dune::GridGlue::GridGlue
sequential adapter to couple two grids at specified close together boundaries
Definition: gridglue.hh:40
Dune::GridGlue::Reverse
Definition: rangegenerators.hh:13