Moertel  Development
 All Classes Namespaces Files Functions Enumerations Friends Pages
Classes | Functions
MOERTEL Namespace Reference

MOERTEL: namespace of the Moertel package. More...

Classes

class  Function
 A virtual class to form a shape function of some type. More...
 
class  Function_Constant1D
  A 1D constant shape function for a 2-noded 1D Segment More...
 
class  Function_Linear1D
  A 1D linear shape function More...
 
class  Function_DualLinear1D
  A 1D dual linear shape function More...
 
class  Function_LinearTri
  A 2D linear shape function for a 3-noded triangle More...
 
class  Function_DualLinearTri
  A 2D dual linear shape function for a 3-noded triangle More...
 
class  Function_ConstantTri
  A 2D constant shape function for a 3-noded triangle More...
 
class  Function_BiLinearQuad
 
class  Function_DualBiLinearQuad
 
class  Integrator
 A class to perform integration of mass matrices on the overlap of 2 segments in 1D and 2D. More...
 
class  Interface
 A class to construct a single interface. More...
 
class  Lmselector
 A virtual class that supports user activation of LMs. More...
 
class  Manager
 Top level user interface to the mortar package 'Moertel'. More...
 
class  Mortar_ML_Preconditioner
 
class  Node
 A class to construct a single node. More...
 
class  ProjectedNode
  A class to handle the projection of a node onto some segment More...
 
class  Point
 A light weight node. More...
 
class  Projector
 A class to perform projections of nodes onto opposing segments in 2D and 3D. More...
 
class  Segment
 A virtual class to serve as base class for different types of interface segmentsconstruct a single interface. More...
 
class  Segment_BiLinearQuad
 A class to define a 4-noded quadrilateral 2D Segment. More...
 
class  Segment_BiLinearTri
 A class to define a 3-noded triangle 2D Segment. More...
 
class  Segment_Linear1D
 A class to define a 2-noded linear 1D Segment. More...
 
class  Solver
 A class to solve mortar constraint problems. More...
 

Functions

MOERTEL::FunctionAllocateFunction (MOERTEL::Function::FunctionType type, int out)
 Allocate a function of the correct type. More...
 
MOERTEL::SegmentAllocateSegment (int type, int out)
 Allocate a Segment of the correct type. More...
 
bool cross (double *out, const double *g1, const double *g2)
 Cross product. More...
 
double dot (const double *g1, const double *g2, const int dim)
 Dot product. More...
 
double length (const double *g, const int dim)
 Length of a vector. More...
 
bool solve22 (const double A[][2], double *x, const double *b)
 Solve dense 2x2 system of equations. More...
 
bool solve33 (const double A[][3], double *x, const double *b)
 Solve dense 3x3 system of equations. More...
 
int digit_ten (int i)
 Return the '10' digit from an integer number.
 
void sort (double *dlist, int N, int *list2)
 Sort dlist of length N in ascending, sort list2 according to dlist. More...
 
template<typename kind >
void swap (kind &a, kind &b)
 Template to swap 2 <kind> instances. More...
 
int MatrixMatrixAdd (const Epetra_CrsMatrix &A, bool transposeA, double scalarA, Epetra_CrsMatrix &B, double scalarB)
 Add matrices A+B. More...
 
Epetra_CrsMatrix * MatMatMult (const Epetra_CrsMatrix &A, bool transA, const Epetra_CrsMatrix &B, bool transB, int outlevel)
 Multiply matrices A*B. More...
 
Epetra_CrsMatrix * PaddedMatrix (const Epetra_Map rowmap, double val, const int numentriesperrow)
 Allocate and return a matrix padded with val on the diagonal. FillComplete() is NOT called on exit.
 
Epetra_CrsMatrix * StripZeros (Epetra_CrsMatrix &A, double eps)
 Strip out values from a matrix below a certain tolerance. More...
 
bool SplitMatrix2x2 (Teuchos::RCP< Epetra_CrsMatrix > A, Teuchos::RCP< Epetra_Map > &A11rowmap, Teuchos::RCP< Epetra_Map > &A22rowmap, Teuchos::RCP< Epetra_CrsMatrix > &A11, Teuchos::RCP< Epetra_CrsMatrix > &A12, Teuchos::RCP< Epetra_CrsMatrix > &A21, Teuchos::RCP< Epetra_CrsMatrix > &A22)
 split a matrix into a 2x2 block system where the rowmap of one of the blocks is given More...
 
Epetra_Map * SplitMap (const Epetra_Map &Amap, const Epetra_Map &Agiven)
 split a rowmap of matrix A More...
 
bool SplitVector (const Epetra_Vector &x, const Epetra_Map &x1map, Epetra_Vector *&x1, const Epetra_Map &x2map, Epetra_Vector *&x2)
 split a vector into 2 non-overlapping pieces
 
bool MergeVector (const Epetra_Vector &x1, const Epetra_Vector &x2, Epetra_Vector &xresult)
 merge results from 2 vectors into one (assumes matching submaps)
 
bool Print_Matrix (std::string name, Epetra_CrsMatrix &A, int ibase)
 Print matrix to file. More...
 
bool Print_Graph (std::string name, Epetra_CrsGraph &A, int ibase)
 Print graph to file. More...
 
bool Print_Vector (std::string name, Epetra_Vector &v, int ibase)
 Print vector to file. More...
 
int ReportError (const std::stringstream &Message)
 Error reporting method.
 

Detailed Description

MOERTEL: namespace of the Moertel package.

The Moertel package depends on Epetra, EpetraExt, Teuchos, Amesos, ML and AztecOO:
Use at least the following lines in the configure of Trilinos:

--enable-moertel
--enable-epetra
--enable-epetraext
--enable-teuchos
--enable-ml
--enable-aztecoo --enable-aztecoo-teuchos
--enable-amesos

Function Documentation

MOERTEL::Function * MOERTEL::AllocateFunction ( MOERTEL::Function::FunctionType  type,
int  out 
)

Allocate a function of the correct type.

For communication reasons, every single derived function class needs to have a unique typ-id. This type Id can be communicated easily. So when introducing a new derived Function class, one needs to add it's type to the enum FunctionType in the virtual base class in mrtr_function.H and one needs to add a case to this method MOERTEL::AllocateFunction in mrtr_utils.cpp

Parameters
type: Type of Function to allocate and return pointer to

References ReportError().

Referenced by MOERTEL::Segment_Linear1D::UnPack(), MOERTEL::Segment_BiLinearTri::UnPack(), and MOERTEL::Segment_BiLinearQuad::UnPack().

MOERTEL::Segment * MOERTEL::AllocateSegment ( int  type,
int  out 
)

Allocate a Segment of the correct type.

For communication reasons, every single derived segment class needs to have a unique typ-id. This type Id can be communicated easily. So when introducing a new segment class, one needs to add it's type to the enum SegmentType in the virtual base class in mrtr_segment.H and one needs to add a case to this method MOERTEL::AllocateSegment in mrtr_utils.cpp

Parameters
type: Type of segment to allocate and return pointer to
out: Level of output to be generated to stdout ( 0 - 10 )

References ReportError().

bool MOERTEL::cross ( double *  out,
const double *  g1,
const double *  g2 
)

Cross product.

Perform cross product out = g1 x g2 for vectors of dimension 3

Referenced by MOERTEL::Segment_BiLinearTri::BuildNormal(), and MOERTEL::Segment_BiLinearQuad::Metric().

double MOERTEL::dot ( const double *  g1,
const double *  g2,
const int  dim 
)

Dot product.

Perform dot product g1 dot g2 for vectors of dimension dim and return result

Referenced by MOERTEL::Projector::ProjectNodetoSegment_NodalNormal(), and MOERTEL::Projector::ProjectNodetoSegment_SegmentNormal().

double MOERTEL::length ( const double *  g,
const int  dim 
)
Epetra_CrsMatrix * MOERTEL::MatMatMult ( const Epetra_CrsMatrix &  A,
bool  transA,
const Epetra_CrsMatrix &  B,
bool  transB,
int  outlevel 
)

Multiply matrices A*B.

matrices A and B are mutliplied and the result is allocated and returned. The method makes uses EpetraExt for multiplication The user is responsible for freeing the returned result.

Parameters
A: Matrix A to multiply
transA: flag indicating whether A*T shall be used
B: Matrix B to multiply
transB: flag indicating whether B*T shall be used
Returns
Result upon success and NULL upon failure

Referenced by MOERTEL::Mortar_ML_Preconditioner::Compute(), and MOERTEL::Manager::MakeSPDProblem().

int MOERTEL::MatrixMatrixAdd ( const Epetra_CrsMatrix &  A,
bool  transposeA,
double  scalarA,
Epetra_CrsMatrix &  B,
double  scalarB 
)

Add matrices A+B.

Perform B = scalarB * B + scalarA * A ^ transposeA If scalarB is 0.0, then B = scalarA * A ^ transposeA isperformed.

This is a modified version of EpetraExt's MatrixMatrixAdd. FillComplete() must not be called on B upon entry, FillComplete() will not be called on B upon exit by this method.

Parameters
A: Matrix A to add to B
transposeA: flag indicating whether A*T shall be added
scalarA: scalar factor for A
B: Matrix B to be added to
scalarB: scalar factor for B
Returns
Zero upon success

References ReportError().

Referenced by MOERTEL::Mortar_ML_Preconditioner::Compute(), MOERTEL::Manager::MakeSaddleProblem(), and MOERTEL::Manager::MakeSPDProblem().

bool MOERTEL::Print_Graph ( std::string  name,
Epetra_CrsGraph &  A,
int  ibase 
)

Print graph to file.

Prints an Epetra_CrsGraph to file in serial and parallel. Will create several files with process id appended to the name in parallel. Index base can either be 0 or 1. The first row of the file gives the global size of the range and domain map, the second row gives the local size of the row- and column map.

Parameters
name: Name of file without appendix, appendix will be .mtx
A: Graph to print
ibase: Index base, should be either 1 or 0

References ReportError().

bool MOERTEL::Print_Matrix ( std::string  name,
Epetra_CrsMatrix &  A,
int  ibase 
)

Print matrix to file.

Prints an Epetra_CrsMatrix to file in serial and parallel. Will create several files with process id appended to the name in parallel. Index base can either be 0 or 1. The first row of the file gives the global size of the range and domain map, the sond row gives the local size of the row- and column map.

Parameters
name: Name of file without appendix, appendix will be .mtx
A: Matrix to print
ibase: Index base, should be either 1 or 0

References ReportError().

bool MOERTEL::Print_Vector ( std::string  name,
Epetra_Vector &  v,
int  ibase 
)

Print vector to file.

Prints an Epetra_Vector to file in serial and parallel. Will create several files with process id appended to the name in parallel. Index base can either be 0 or 1.

Parameters
name: Name of file without appendix, appendix will be .vec
v: Vector to print
ibase: Index base, should be either 1 or 0
bool MOERTEL::solve22 ( const double  A[][2],
double *  x,
const double *  b 
)

Solve dense 2x2 system of equations.

Ax=b

References ReportError().

bool MOERTEL::solve33 ( const double  A[][3],
double *  x,
const double *  b 
)
void MOERTEL::sort ( double *  dlist,
int  N,
int *  list2 
)

Sort dlist of length N in ascending, sort list2 according to dlist.

This piece of code was lend from the Trilinos package ML

Epetra_Map * MOERTEL::SplitMap ( const Epetra_Map &  Amap,
const Epetra_Map &  Agiven 
)

split a rowmap of matrix A

splits A->RowMap() into 2 maps and returns them, where one of the rowmaps has to be given on input

Parameters
Amap: Map to split on input
Agiven: on entry submap that is given and part of Amap
Returns
the remainder map of Amap that is not overlapping with Agiven

Referenced by SplitMatrix2x2().

bool MOERTEL::SplitMatrix2x2 ( Teuchos::RCP< Epetra_CrsMatrix >  A,
Teuchos::RCP< Epetra_Map > &  A11rowmap,
Teuchos::RCP< Epetra_Map > &  A22rowmap,
Teuchos::RCP< Epetra_CrsMatrix > &  A11,
Teuchos::RCP< Epetra_CrsMatrix > &  A12,
Teuchos::RCP< Epetra_CrsMatrix > &  A21,
Teuchos::RCP< Epetra_CrsMatrix > &  A22 
)

split a matrix into a 2x2 block system where the rowmap of one of the blocks is given

Splits a given matrix into a 2x2 block system where the rowmap of one of the blocks is given on input. Blocks A11 and A22 are assumed to be square. All values on entry have to be Teuchos::null except the given rowmap and matrix A. Note that either A11rowmap or A22rowmap or both have to be nonzero. In case both rowmaps are supplied they have to be an exact and nonoverlapping split of A->RowMap(). Matrix blocks are FillComplete() on exit.

Parameters
A: Matrix A on input
A11rowmap: rowmap of A11 or null
A22rowmap: rowmap of A22 or null
A11: on exit matrix block A11
A12: on exit matrix block A12
A21: on exit matrix block A21
A22: on exit matrix block A22

References length(), ReportError(), and SplitMap().

Epetra_CrsMatrix * MOERTEL::StripZeros ( Epetra_CrsMatrix &  A,
double  eps 
)

Strip out values from a matrix below a certain tolerance.

Allocates and returns a new matrix and copies A to it where entries with an absoute value smaller then eps are negelected. The method calls FillComplete(A.OperatorDomainMap(),A.OperatorRangeMap()) on the result.

Parameters
A: Matrix A to strip
eps: tolerance
Returns
The new matrix upon success, NULL otherwise

References ReportError().

template<typename kind >
void MOERTEL::swap ( kind &  a,
kind &  b 
)

Template to swap 2 <kind> instances.

<kind> has to implement the assignment operator =