Top level user interface to the mortar package 'Moertel'. More...
#include <mrtr_manager.H>
Public Types | |
enum | DimensionType { manager_none, manager_2D, manager_3D } |
Specification of the problem dimension. | |
Public Member Functions | |
Manager (Epetra_Comm &comm, MOERTEL::Manager::DimensionType dimension, int outlevel) | |
Creates and empty instance of this class. More... | |
Manager (Epetra_Comm &comm, int outlevel) | |
Creates and empty instance of this class. More... | |
virtual | ~Manager () |
Destroys an instance of this class. More... | |
Epetra_Comm & | Comm () const |
Returns the Epetra_Comm object associated with this class. | |
int | OutLevel () |
Returns the Level of output (0 - 10) the user specified in the constructor. | |
MOERTEL::Manager::DimensionType | Dimension () |
Query the problem dimension. More... | |
int | Ninterfaces () |
Query the number of interfaces passed to this class via AddInterface. More... | |
bool | Print () const |
Print all information stored in this class to stdout. | |
void | SetDimension (MOERTEL::Manager::DimensionType type) |
Set problem dimension. More... | |
bool | AddInterface (MOERTEL::Interface &interface) |
Add an interface class to this class. More... | |
Teuchos::ParameterList & | Default_Parameters () |
Obtain a parameter list with default values. More... | |
bool | Mortar_Integrate () |
Perform integration of the mortar integral on all interfaces. More... | |
bool | Integrate_Interfaces () |
Perform the integration of the mortar integral on all interfaces. More... | |
bool | AssembleInterfacesIntoResidual (Lmselector *sel) |
Assembles contributions from the D and M matrices into the JFNK residual vector. More... | |
bool | SetProblemMap (const Epetra_Map *map) |
Set the RowMap of the original uncoupled problem. More... | |
Epetra_Map * | ProblemMap () const |
Return view of row map of the uncoupled problem. More... | |
bool | SetInputMatrix (Epetra_CrsMatrix *inputmatrix, bool DeepCopy=false) |
Set the Matrix of the original uncoupled problem. More... | |
Epetra_CrsMatrix * | MakeSaddleProblem () |
Construct a saddle point system of equations. More... | |
Epetra_CrsMatrix * | MakeSPDProblem () |
Construct a spd system of equations if possible. More... | |
Epetra_CrsMatrix * | GetRHSMatrix () |
returns the righ hand side matrix of the spd system of equations More... | |
void | SetSolverParameters (Teuchos::ParameterList ¶ms) |
Set a parameter list holding solver parameters. More... | |
bool | Solve (Teuchos::ParameterList ¶ms, Epetra_Vector &sol, const Epetra_Vector &rhs) |
Solve the problem. More... | |
bool | Solve (Epetra_Vector &sol, const Epetra_Vector &rhs) |
Solve the problem. More... | |
void | ResetSolver () |
Reset the internal solver. More... | |
const Epetra_CrsMatrix * | D () const |
Get pointer to constraints matrix D. | |
const Epetra_CrsMatrix * | M () const |
Get pointer to constraints matrix M. | |
std::map< int, int > | lm_to_dof () |
Teuchos::RCP< MOERTEL::Interface > | GetInterface (int idx) |
Query a managed interface added using AddInterface. More... | |
Teuchos::RCP< Epetra_Map > | GetSaddleMap () |
Protected Attributes | |
int | outlevel_ |
Epetra_Comm & | comm_ |
DimensionType | dimensiontype_ |
std::map< int, Teuchos::RCP < MOERTEL::Interface > > | interface_ |
Teuchos::RCP< Epetra_Map > | problemmap_ |
Teuchos::RCP< Epetra_CrsMatrix > | inputmatrix_ |
Teuchos::RCP< Epetra_Map > | constraintsmap_ |
Teuchos::RCP< Epetra_CrsMatrix > | D_ |
Teuchos::RCP< Epetra_CrsMatrix > | M_ |
Teuchos::RCP< Epetra_Map > | saddlemap_ |
Teuchos::RCP< Epetra_CrsMatrix > | saddlematrix_ |
Teuchos::RCP< Epetra_CrsMatrix > | spdmatrix_ |
Teuchos::RCP< Epetra_CrsMatrix > | spdrhs_ |
Teuchos::RCP< std::map< int, int > > | lm_to_dof_ |
Teuchos::RCP< Epetra_CrsMatrix > | I_ |
Teuchos::RCP< Epetra_CrsMatrix > | WT_ |
Teuchos::RCP< Epetra_CrsMatrix > | B_ |
Teuchos::RCP< Epetra_Map > | annmap_ |
Teuchos::RCP < Teuchos::ParameterList > | integrationparams_ |
Teuchos::RCP < Teuchos::ParameterList > | solverparams_ |
Teuchos::RCP< MOERTEL::Solver > | solver_ |
Friends | |
class | MOERTEL::Solver |
Top level user interface to the mortar package 'Moertel'.
Top level user interface to mortar package 'Moertel'
This class supplies capabilities for nonconforming mesh tying in 2 and 3 dimensions using Mortar methods.
It also supports use of mortar methods for contact formulations and domain decomposition as well as some limited support to solving mortar-coupled systems of equations.
Proceeding:
The MOERTEL::Manager class supports the std::ostream& operator <<
The package and this class make use of the other Trilinos packages Epetra , EpetraExt , Teuchos , Amesos , ML .
For more detailed information on the usage, see "The Moertel User's Guide", Sandia report XXXX
For background information on Mortar methods, we refer to
Wohlmuth, B.:"Discretization Methods and Iterative Solvers Based on Domain Decomposition", Springer, New York, 2001, ISBN 3-540-41083-X
|
explicit |
Creates and empty instance of this class.
Constructs an empty instance of this class to be filled by the user with information about the non-conforming interface(s).
This is a collective call for all processors associated with the Epetra_Comm
comm | : An Epetra_Comm object |
dimension | : Dimension of the problem |
outlevel | : The level of output to stdout to be generated ( 0 - 10 ) |
References Default_Parameters().
|
explicit |
Creates and empty instance of this class.
Constructs an empty instance of this class that is filled in by the user with information about the non-conforming interface(s).
This is a collective call for all processors associated with the Epetra_Comm
comm | : An Epetra_Comm object |
outlevel | : The level of output to stdout to be generated ( 0 - 10 ) |
References Default_Parameters().
|
virtual |
Destroys an instance of this class.
Destructor
bool MOERTEL::Manager::AddInterface | ( | MOERTEL::Interface & | interface | ) |
Add an interface class to this class.
Add a previously constructed interface to this class. Before adding an interface, the interface's public method MOERTEL::Interface::Complete() has to be called. This class will not accept interface that are not completed. Any number of interfaces can be added. This class does not take ownership over the added class. The interface added can be destroyed immediately after the call to this method as this class stores a deep copy of it.
interface | : A completed interface |
References MOERTEL::Interface::IsComplete().
bool MOERTEL::Manager::AssembleInterfacesIntoResidual | ( | Lmselector * | sel | ) |
Assembles contributions from the D and M matrices into the JFNK residual vector.
Once all interfaces and a row map of the original system are passed to this class this method evaluates the mortar integrals on all interfaces.
sel | (in): A user-derived object inheriting from Lmselector. |
rhs | (out) : JFNK rhs vector containing accumulated D and M values |
soln | (in) : Current JFNK solution vector, used to matvec with D and M |
Teuchos::ParameterList & MOERTEL::Manager::Default_Parameters | ( | ) |
Obtain a parameter list with default values.
Obtain a parameter list containing default values for integration of mortar interfaces. Default values set are:
Referenced by Manager().
|
inline |
Query the problem dimension.
Query the problem dimension
|
inline |
Query a managed interface added using AddInterface.
Accessor function to the interfaces stored within the manager by previous AddInterface calls.
idx | : index of the interface to be returned (0 up to Ninterfaces) |
Epetra_CrsMatrix * MOERTEL::Manager::GetRHSMatrix | ( | ) |
returns the righ hand side matrix of the spd system of equations
If dual Lagrange multipliers are used, a symm. positive definite system of equations can be constructed from the saddle point problem. The final solution of the problem is then obtained from
Atilde x = RHSMatrix * b
bool MOERTEL::Manager::Integrate_Interfaces | ( | ) |
Perform the integration of the mortar integral on all interfaces.
Once all interfaces and a row map of the original system are passed to this class this method evaluates the mortar integrals on all interfaces.
Process:
Epetra_CrsMatrix * MOERTEL::Manager::MakeSaddleProblem | ( | ) |
Construct a saddle point system of equations.
After a call to Mortar_Integrate() a saddle point system of equations can be constructed and returned to the user. Prerequisite is that the user has supplied the original uncoupled matrix via SetInputMatrix and the integration was performed with Mortar_Integrate(). This class holds ownership of the saddle point system so the user must not destroy it.
References MOERTEL::Solver::Comm(), MOERTEL::MatrixMatrixAdd(), and MOERTEL::Solver::OutLevel().
Epetra_CrsMatrix * MOERTEL::Manager::MakeSPDProblem | ( | ) |
Construct a spd system of equations if possible.
If dual Lagrange multipliers are used, a symm. positive definite system of equations can be constructed from the saddle point problem. The final solution of the problem is then obtained from
Atilde x = RHSMatrix * b
References MOERTEL::Solver::Comm(), MOERTEL::MatMatMult(), MOERTEL::MatrixMatrixAdd(), MOERTEL::Solver::OutLevel(), and MOERTEL::PaddedMatrix().
bool MOERTEL::Manager::Mortar_Integrate | ( | ) |
Perform integration of the mortar integral on all interfaces.
Once all interfaces and a row map of the original system are passed to this class this method evaluates the mortar integrals on all interfaces.
Process:
|
inline |
Query the number of interfaces passed to this class via AddInterface.
Returns the number of non-conforming interfaces passed to this class via AddInterface
|
inline |
Return view of row map of the uncoupled problem.
return a view of the row map of the uncoupled problem the user provided using SetProblemMap. Returns NULL if no row map was previously set.
void MOERTEL::Manager::ResetSolver | ( | ) |
Reset the internal solver.
The internal solver (no matter which one) will go on using the once constructed (and maybe factored) matrix until this method is called. After a call to this method, the linear problem as well as the solver are build from scratch.
|
inline |
Set problem dimension.
This class can handle 2D and 3D problems but not both at the same time. It is necessary to specify the dimension of the problem
type | : Dimension of the problem |
bool MOERTEL::Manager::SetInputMatrix | ( | Epetra_CrsMatrix * | inputmatrix, |
bool | DeepCopy = false |
||
) |
Set the Matrix of the original uncoupled problem.
Passes in the Matrix of the original (uncoupled) problem. In a saddle point system, this would be the (1,1) block. This class takes ownership of the matrix if DeepCopy is set to true.
The matrix of the original, uncoupled problem is used to construct a saddle point system via MakeSaddleProblem() or solve the coupled problem
map | : Matrix of the original (uncoupled) problem |
bool MOERTEL::Manager::SetProblemMap | ( | const Epetra_Map * | map | ) |
Set the RowMap of the original uncoupled problem.
Passes in the RowMap of the original (uncoupled) problem. In a saddle point system, this would be the row map of the (1,1) block. This class does not take ownership of the passed in map. The map can be destroyed immediately after the call to this method as this class stores a deep copy of it. A row map has to be passed to this class before a call to Mortar_Integrate()
map | : Row map of the original (uncoupled) problem |
void MOERTEL::Manager::SetSolverParameters | ( | Teuchos::ParameterList & | params | ) |
Set a parameter list holding solver parameters.
Set a ParameterList containing solver options. This class does not take ownership over params but instead uses a view of it.
Currently Moertel has interfaces to the Amesos, the AztecOO and the ML package.
All packages are controlled by this parameter list. The parameter list Contains a general part where the package and the type of linear system to be generated is chosen and contains one or more sub-parameter lists which are passed on to the appropriate package. This way all parameters described in the documentation of Amesos, ML and Aztec can be passed to the appropriate package.
Below one example how to choose the parameter list is given:
params | : Parameter list containing solver options |
bool MOERTEL::Manager::Solve | ( | Teuchos::ParameterList & | params, |
Epetra_Vector & | sol, | ||
const Epetra_Vector & | rhs | ||
) |
Solve the problem.
Solve the coupled problem using the solver options supplied in params. Succesive calls to this method will reuse the system matrix together with the supplied solution and rhs vectors until ResetSolver() is called. Then, the linear problem is recreated from scratch.
params | (in) : Solution parameters |
sol | (out) : The solution |
rhs | (in) : The right hand side vector |
bool MOERTEL::Manager::Solve | ( | Epetra_Vector & | sol, |
const Epetra_Vector & | rhs | ||
) |
Solve the problem.
Solve the coupled problem using the solver options previously supplied in using SetSolverParameters Succesive calls to this method will resue the system matrix together with the supplied solution and rhs vectors until ResetSolver() is called. Then, the linear problem is recreated from scratch.
sol | (out) : The solution |
rhs | (in) : The right hand side vector |
References MOERTEL::Solver::Comm(), and MOERTEL::Solver::OutLevel().