Moertel  Development
 All Classes Namespaces Files Functions Enumerations Friends Pages
List of all members
MOERTEL::Integrator Class Reference

A class to perform integration of mass matrices on the overlap of 2 segments in 1D and 2D. More...

#include <mrtr_integrator.H>

Public Member Functions

 Integrator (int ngp, bool oneD, int outlevel)
 Constructor. More...
 
virtual ~Integrator ()
 Destructor.
 
int OutLevel ()
 Return the level of output written to stdout ( 0 - 10 )
 
int Ngp ()
 Return number of Gaussian integration points to be used.
 
double Coordinate (int gp)
 Return coordinate of a specific Gaussian point in 1D segment. More...
 
double * Coordinate (int *gp)
 Return coordinates of a specific Gaussian point in 2D segment. More...
 
double Weight (int gp)
 Return weight for given Gaussian point. More...
 
Epetra_SerialDenseMatrix * Integrate (MOERTEL::Segment &sseg, double sxia, double sxib, MOERTEL::Segment &mseg, double mxia, double mxib)
 Integrate mass matrix 'M' on a 1D overlap between a slave and a master segment. More...
 
bool Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_CrsMatrix &M, Epetra_SerialDenseMatrix &Mdense)
 Assemble integration result '-M' into global matrix 'M'. More...
 
Epetra_SerialDenseMatrix * Integrate (MOERTEL::Segment &sseg, double sxia, double sxib)
 Integrate mass matrix 'D' on a 1D slave segment overlap. More...
 
bool Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, Epetra_CrsMatrix &D, Epetra_SerialDenseMatrix &Ddense)
 Assemble integration result 'D' into global matrix 'D'. More...
 
Epetra_SerialDenseMatrix * Integrate_2D_Mmod (MOERTEL::Segment &sseg, double sxia, double sxib, MOERTEL::Segment &mseg, double mxia, double mxib)
 Integrate modification to mass matrix 'M' on a 1D overlap between a slave and a master segment. More...
 
bool Assemble_2D_Mod (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_SerialDenseMatrix &Mmod)
 Assemble modification integration result '-DELTA_M' into global matrix 'M'. More...
 
bool Integrate (Teuchos::RCP< MOERTEL::Segment > actseg, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_SerialDenseMatrix **Ddense, Epetra_SerialDenseMatrix **Mdense, MOERTEL::Overlap< MOERTEL::Interface > &overlap, double eps, bool exactvalues)
 Integrate a 2D triangle overlap segment, master part M and slave part D. More...
 
bool Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, Epetra_SerialDenseMatrix &Ddense)
 Assemble integration result 'D' into Node (2D interfaces only) More...
 
bool Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_SerialDenseMatrix &Mdense)
 Assemble integration result 'M' into Node (2D interfaces only) More...
 

Detailed Description

A class to perform integration of mass matrices on the overlap of 2 segments in 1D and 2D.

A class to perform Gaussian integration and assembly of mass matrices on the overlap of 2 segments in 1D and 2D

Date
Last update do Doxygen: 20-March-06
Author
Glen Hansen (gahan.nosp@m.se@s.nosp@m.andia.nosp@m..gov)

Constructor & Destructor Documentation

MOERTEL::Integrator::Integrator ( int  ngp,
bool  oneD,
int  outlevel 
)
explicit

Constructor.

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

Parameters
ngp: number of gaussian points to be used
oneD: flag indicating whether 1D or 2D regions shall be integrated
outlevel: Level of output information written to stdout ( 0 - 10 )

References MOERTEL::ReportError().

Member Function Documentation

bool MOERTEL::Integrator::Assemble ( MOERTEL::Interface inter,
MOERTEL::Segment sseg,
MOERTEL::Segment mseg,
Epetra_CrsMatrix &  M,
Epetra_SerialDenseMatrix &  Mdense 
)

Assemble integration result '-M' into global matrix 'M'.

Parameters
inter: Interface
sseg: Slave Segment
mseg: Mortar Segment
M: global sparse matrix 'M'
Mdense: local dense matrix from integration of overlap between sseg and mseg

References MOERTEL::Node::Dof(), MOERTEL::Node::Id(), MOERTEL::Interface::lComm(), MOERTEL::Node::LMDof(), MOERTEL::Node::Ndof(), MOERTEL::Node::Nlmdof(), MOERTEL::Segment::Nnode(), MOERTEL::Interface::NodePID(), MOERTEL::Segment::Nodes(), and MOERTEL::ReportError().

bool MOERTEL::Integrator::Assemble ( MOERTEL::Interface inter,
MOERTEL::Segment sseg,
Epetra_CrsMatrix &  D,
Epetra_SerialDenseMatrix &  Ddense 
)

Assemble integration result 'D' into global matrix 'D'.

Parameters
inter: Interface
sseg: Slave Segment
D: global sparse matrix 'D'
Ddense: local dense matrix from integration of overlap between sseg and mseg

References MOERTEL::Node::Dof(), MOERTEL::Node::Id(), MOERTEL::Interface::lComm(), MOERTEL::Node::LMDof(), MOERTEL::Node::Ndof(), MOERTEL::Node::Nlmdof(), MOERTEL::Segment::Nnode(), MOERTEL::Interface::NodePID(), MOERTEL::Segment::Nodes(), and MOERTEL::ReportError().

bool MOERTEL::Integrator::Assemble ( MOERTEL::Interface inter,
MOERTEL::Segment sseg,
Epetra_SerialDenseMatrix &  Ddense 
)

Assemble integration result 'D' into Node (2D interfaces only)

Parameters
inter: Interface
sseg: Slave Segment
Ddense: local dense matrix from integration of overlap

References MOERTEL::Node::AddMValue(), MOERTEL::Node::GetSupportedByNode(), MOERTEL::Node::Id(), MOERTEL::Interface::lComm(), MOERTEL::Segment::Nnode(), MOERTEL::Interface::NodePID(), and MOERTEL::Segment::Nodes().

bool MOERTEL::Integrator::Assemble ( MOERTEL::Interface inter,
MOERTEL::Segment sseg,
MOERTEL::Segment mseg,
Epetra_SerialDenseMatrix &  Mdense 
)

Assemble integration result 'M' into Node (2D interfaces only)

Parameters
inter: Interface
sseg: Slave Segment
mseg: Mortar Segment
Mdense: local dense matrix from integration of overlap

References MOERTEL::Node::GetSupportedByNode(), MOERTEL::Node::Id(), MOERTEL::Interface::lComm(), MOERTEL::Segment::Nnode(), MOERTEL::Interface::NodePID(), and MOERTEL::Segment::Nodes().

bool MOERTEL::Integrator::Assemble_2D_Mod ( MOERTEL::Interface inter,
MOERTEL::Segment sseg,
MOERTEL::Segment mseg,
Epetra_SerialDenseMatrix &  Mmod 
)

Assemble modification integration result '-DELTA_M' into global matrix 'M'.

Parameters
inter: Interface
sseg: Slave Segment
mseg: Mortar Segment
M: global sparse matrix 'M'
Mmod: local dense matrix from integration of modification on overlap between sseg and mseg

References MOERTEL::Node::Dof(), MOERTEL::Node::GetSupportedByNode(), MOERTEL::Node::Id(), MOERTEL::Interface::lComm(), MOERTEL::Node::Ndof(), MOERTEL::Segment::Nnode(), MOERTEL::Interface::NodePID(), and MOERTEL::Segment::Nodes().

double MOERTEL::Integrator::Coordinate ( int  gp)
inline

Return coordinate of a specific Gaussian point in 1D segment.

Parameters
gp: Number of Gaussian point to get the coordinate for
double* MOERTEL::Integrator::Coordinate ( int *  gp)
inline

Return coordinates of a specific Gaussian point in 2D segment.

Parameters
gp: Number of Gaussian point to get the coordinate for
Epetra_SerialDenseMatrix * MOERTEL::Integrator::Integrate ( MOERTEL::Segment sseg,
double  sxia,
double  sxib,
MOERTEL::Segment mseg,
double  mxia,
double  mxib 
)

Integrate mass matrix 'M' on a 1D overlap between a slave and a master segment.

Integrate over overlap the trace space shape function of the mortar side times the Lagrange multiplier space function of the slave side

Parameters
sseg: Slave Segment
sxia: lower bound of overlap in slave segment coordinates
sxib: upper bound of overlap in slave segment coordinates
mseg: Mortar Segment
mxia: lower bound of overlap in mortar segment coordinates
mxib: upper bound of overlap in mortar segment coordinates

References MOERTEL::Segment::EvaluateFunction(), MOERTEL::Segment::Metric(), and MOERTEL::Segment::Nnode().

Epetra_SerialDenseMatrix * MOERTEL::Integrator::Integrate ( MOERTEL::Segment sseg,
double  sxia,
double  sxib 
)

Integrate mass matrix 'D' on a 1D slave segment overlap.

Integrate over overlap the trace space shape function of the slave side times the Lagrange multiplier space function of the slave side

Parameters
sseg: Slave Segment
sxia: lower bound of overlap in slave segment coordinates
sxib: upper bound of overlap in slave segment coordinates

References MOERTEL::Segment::EvaluateFunction(), MOERTEL::Segment::Metric(), and MOERTEL::Segment::Nnode().

bool MOERTEL::Integrator::Integrate ( Teuchos::RCP< MOERTEL::Segment actseg,
MOERTEL::Segment sseg,
MOERTEL::Segment mseg,
Epetra_SerialDenseMatrix **  Ddense,
Epetra_SerialDenseMatrix **  Mdense,
MOERTEL::Overlap< MOERTEL::Interface > &  overlap,
double  eps,
bool  exactvalues 
)

Integrate a 2D triangle overlap segment, master part M and slave part D.

Integrate a triangle which is part of the discretization of a overlap polygon between a slave and a mortar segment.

Parameters
actseg: triangle to be integrated over
sseg: Slave Segment
mseg: Mortar Segment
Ddense: Dense matrix holding integration result 'D'
Mdense: Dense matrix holding integration result 'M'
overlap: the overlap class holds all neccessary function values
eps: Threshold, skip triangles with an area < eps

References MOERTEL::Segment::Area(), MOERTEL::Segment::EvaluateFunction(), MOERTEL::length(), MOERTEL::Segment::Nnode(), MOERTEL::Projector::ProjectNodetoSegment_NodalNormal(), and MOERTEL::ReportError().

Epetra_SerialDenseMatrix * MOERTEL::Integrator::Integrate_2D_Mmod ( MOERTEL::Segment sseg,
double  sxia,
double  sxib,
MOERTEL::Segment mseg,
double  mxia,
double  mxib 
)

Integrate modification to mass matrix 'M' on a 1D overlap between a slave and a master segment.

Integrate over overlap the modification of the trace space shape function of the mortar side times the Lagrange multiplier space function of the slave side. This modification due to B. wohlmuth improves behaviour for curved interfaces

Parameters
sseg: Slave Segment
sxia: lower bound of overlap in slave segment coordinates
sxib: upper bound of overlap in slave segment coordinates
mseg: Mortar Segment
mxia: lower bound of overlap in mortar segment coordinates
mxib: upper bound of overlap in mortar segment coordinates

References MOERTEL::Segment::EvaluateFunction(), MOERTEL::Segment::Metric(), and MOERTEL::Segment::Nnode().

double MOERTEL::Integrator::Weight ( int  gp)
inline

Return weight for given Gaussian point.

Parameters
gp: Number of Gaussian point to get the coordinate for

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