#include <mrtr_integrator.H>
Public Member Functions | |
| Integrator (int ngp, bool oneD, int outlevel) | |
| Constructor. | |
| 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. | |
| double * | Coordinate (int *gp) |
| Return coordinates of a specific Gaussian point in 2D segment. | |
| double | Weight (int gp) |
| Return weight for given Gaussian point. | |
| 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. | |
| 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'. | |
| Epetra_SerialDenseMatrix * | Integrate (MOERTEL::Segment &sseg, double sxia, double sxib) |
| Integrate mass matrix 'D' on a 1D slave segment overlap. | |
| bool | Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, Epetra_CrsMatrix &D, Epetra_SerialDenseMatrix &Ddense) |
| Assemble integration result 'D' into global matrix 'D'. | |
| 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. | |
| 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'. | |
| bool | Integrate (RefCountPtr< MOERTEL::Segment > actseg, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_SerialDenseMatrix **Ddense, Epetra_SerialDenseMatrix **Mdense, MOERTEL::Overlap &overlap, double eps, bool exactvalues) |
| Integrate a 2D triangle overlap segment, master part M and slave part D. | |
| bool | Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, Epetra_SerialDenseMatrix &Ddense) |
| Assemble integration result 'D' into Node (2D interfaces only). | |
| bool | Assemble (MOERTEL::Interface &inter, MOERTEL::Segment &sseg, MOERTEL::Segment &mseg, Epetra_SerialDenseMatrix &Mdense) |
| Assemble integration result 'M' into Node (2D interfaces only). | |
| 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.
| 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 ) |
| bool MOERTEL::Integrator::Assemble | ( | MOERTEL::Interface & | inter, | |
| MOERTEL::Segment & | sseg, | |||
| MOERTEL::Segment & | mseg, | |||
| Epetra_SerialDenseMatrix & | Mdense | |||
| ) |
| bool MOERTEL::Integrator::Assemble | ( | MOERTEL::Interface & | inter, | |
| MOERTEL::Segment & | sseg, | |||
| Epetra_SerialDenseMatrix & | Ddense | |||
| ) |
| bool MOERTEL::Integrator::Assemble | ( | MOERTEL::Interface & | inter, | |
| MOERTEL::Segment & | sseg, | |||
| Epetra_CrsMatrix & | D, | |||
| Epetra_SerialDenseMatrix & | Ddense | |||
| ) |
| bool MOERTEL::Integrator::Assemble | ( | MOERTEL::Interface & | inter, | |
| MOERTEL::Segment & | sseg, | |||
| MOERTEL::Segment & | mseg, | |||
| Epetra_CrsMatrix & | M, | |||
| Epetra_SerialDenseMatrix & | Mdense | |||
| ) |
| bool MOERTEL::Integrator::Assemble_2D_Mod | ( | MOERTEL::Interface & | inter, | |
| MOERTEL::Segment & | sseg, | |||
| MOERTEL::Segment & | mseg, | |||
| Epetra_SerialDenseMatrix & | Mmod | |||
| ) |
| double* MOERTEL::Integrator::Coordinate | ( | int * | gp | ) | [inline] |
Return coordinates of a specific Gaussian point in 2D segment.
| gp | : Number of Gaussian point to get the coordinate for |
| double MOERTEL::Integrator::Coordinate | ( | int | gp | ) | [inline] |
Return coordinate of a specific Gaussian point in 1D segment.
| gp | : Number of Gaussian point to get the coordinate for |
| bool MOERTEL::Integrator::Integrate | ( | RefCountPtr< MOERTEL::Segment > | actseg, | |
| MOERTEL::Segment & | sseg, | |||
| MOERTEL::Segment & | mseg, | |||
| Epetra_SerialDenseMatrix ** | Ddense, | |||
| Epetra_SerialDenseMatrix ** | Mdense, | |||
| MOERTEL::Overlap & | 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.
| 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 |
| 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
| sseg | : Slave Segment | |
| sxia | : lower bound of overlap in slave segment coordinates | |
| sxib | : upper bound of overlap in slave segment coordinates |
| 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
| 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
| double MOERTEL::Integrator::Weight | ( | int | gp | ) | [inline] |
Return weight for given Gaussian point.
| gp | : Number of Gaussian point to get the coordinate for |
1.4.7