Moertel  Development
Friends | List of all members
MOERTEL::Overlap Class Reference

A class to compute the overlap polygon of 2 different 2D segments and construct a triangle discretization of the convex hull of that polygon. More...

#include <mrtr_overlap.H>

Public Member Functions

 Overlap (MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, MOERTEL::Interface &inter, bool exactvalues, int outlevel)
 Constructor. More...
 
virtual ~Overlap ()
 Destructor.
 
bool ComputeOverlap ()
 Compute overlap (if any) between 2 segments and construct discretization of overlap region. More...
 
int OutLevel ()
 Return the level of output written to stdout ( 0 - 10 )
 

Protected Member Functions

bool AddSegment (int id, MOERTEL::Segment *seg)
 
int Nseg ()
 
void SegmentView (std::vector< Teuchos::RCP< MOERTEL::Segment > > &segs)
 
bool AddPointtoPolygon (const int id, const double *P)
 
bool AddPointtoPolygon (std::map< int, Teuchos::RCP< MOERTEL::Point > > &p, const int id, const double *P)
 
bool RemovePointfromPolygon (const int id, const double *P)
 
int SizePointPolygon ()
 
void PointView (std::vector< Teuchos::RCP< MOERTEL::Point > > &points)
 
void PointView (std::map< int, Teuchos::RCP< MOERTEL::Point > > &p, std::vector< Teuchos::RCP< MOERTEL::Point > > &points)
 
void PointView (std::vector< MOERTEL::Point *> &p, const int *nodeids, const int np)
 
bool CopyPointPolygon (std::map< int, Teuchos::RCP< MOERTEL::Point > > &from, std::map< int, Teuchos::RCP< MOERTEL::Point > > &to)
 
bool Centroid (double xi[], const std::vector< Teuchos::RCP< MOERTEL::Point > > &points, const int np)
 

Friends

class Interface
 the Interface class is a friend to this class
 
class Integrator
 the Integrator class is a friend to this class
 

Detailed Description

A class to compute the overlap polygon of 2 different 2D segments and construct a triangle discretization of the convex hull of that polygon.

A class to compute the overlap polygon of 2 different 2D segments and construct a triangle discretization of the convex hull of that polygon

Date
Last update do Doxygen: 20-March-06

Given a slave and a mortar side segment, this class projects the mortar segment onto the local coordinate system of the slave segment. It will then determine whether an overlap between the slave segment and the projection of the mortar segment exists and computes the polygonal overlap region. In a second step it creates a triangulation of that polygonal overlap region that can be used to perform the integration over that region.

The Interface and the Integrator class are friends to this class to be able to access the resulting triangulation of the overlap polygon.

Author
Glen Hansen (gahan.nosp@m.se@s.nosp@m.andia.nosp@m..gov)

Constructor & Destructor Documentation

◆ Overlap()

MOERTEL::Overlap::Overlap ( MOERTEL::Segment sseg,
MOERTEL::Segment mseg,
MOERTEL::Interface inter,
bool  exactvalues,
int  outlevel 
)
explicit

Constructor.

Constructs an instance of this class.
Note that this is not a collective call as overlaps are computed in parallel by individual processes.

Parameters
sseg: slave Segment to overlap with
mseg: mortar Segment to overlap with
inter: Interface both segments are part of
outlevel: Level of output information written to stdout ( 0 - 10 )

References MOERTEL::ReportError(), and MOERTEL::Segment::Type().

Member Function Documentation

◆ ComputeOverlap()

bool MOERTEL::Overlap::ComputeOverlap ( )

Compute overlap (if any) between 2 segments and construct discretization of overlap region.

Returns
true if segments overlap, false otherwise

The documentation for this class was generated from the following files: