Moertel  Development
 All Classes Namespaces Files Functions Enumerations Friends Pages
Public Types | Protected Attributes | Friends | List of all members
MOERTEL::Manager Class Reference

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 &params)
 Set a parameter list holding solver parameters. More...
 
bool Solve (Teuchos::ParameterList &params, 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::InterfaceGetInterface (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::Solversolver_
 

Friends

class MOERTEL::Solver
 

Detailed Description

Top level user interface to the mortar package 'Moertel'.

Top level user interface to mortar package 'Moertel'

Date
Last update do Doxygen: 09-July-2010

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

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

Constructor & Destructor Documentation

MOERTEL::Manager::Manager ( Epetra_Comm &  comm,
MOERTEL::Manager::DimensionType  dimension,
int  outlevel 
)
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

Parameters
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().

MOERTEL::Manager::Manager ( Epetra_Comm &  comm,
int  outlevel 
)
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

Parameters
comm: An Epetra_Comm object
outlevel: The level of output to stdout to be generated ( 0 - 10 )

References Default_Parameters().

MOERTEL::Manager::~Manager ( )
virtual

Destroys an instance of this class.

Destructor

Member Function Documentation

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.

Parameters
interface: A completed interface
Returns
true when adding the interface was successful and false otherwise

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.

Parameters
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
Returns
true if successful, false otherwise
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:

// When true, values of shape functions are exact at
// every single integration point.
// When false, linear shape functions of an overlap
// discretization are used to linearly interpolate
// shape functions in the overlap region.
// The latter is significantly cheaper but only
// recommended for truly linear functions, e.g. with
// 1D interfaces discretized with linear elements
// 2D interfaces discretized with linear triangle elements
// It is NOT recommended for 2D interfaces discretized with
// bilinear quads as shape functions are not truly linear.
integrationparams_->set("exact values at gauss points",true);
// 1D interface possible values are 1,2,3,4,5,6,7,8,10
integrationparams_->set("number Gaussian points 1D",3);
// 2D interface possible values are 3,6,12,13,16,19,27
integrationparams_->set("number Gaussian points 2D",12);
See Also
SetProblemMap , AddInterface, Mortar_integrate

Referenced by Manager().

MOERTEL::Manager::DimensionType MOERTEL::Manager::Dimension ( )
inline

Query the problem dimension.

Query the problem dimension

Teuchos::RCP<MOERTEL::Interface> MOERTEL::Manager::GetInterface ( int  idx)
inline

Query a managed interface added using AddInterface.

Accessor function to the interfaces stored within the manager by previous AddInterface calls.

Parameters
idx: index of the interface to be returned (0 up to Ninterfaces)
Returns
the requested interface
See Also
AddInterface, 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

See Also
MakeSPDProblem()
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:

  • Chooses a slave and a mortar side for each interface of not already chosen by the user in the interface class
  • Build a C0-continuous field of normals on the slave side
  • Makes special boundary modifications to the Lagrange multiplier shape functions
  • Performs integration of the mortar integrals on individual interfaces
  • Chooses degrees of freedom for the Lagrange multipliers based on the row map of the original uncoupled problem.
  • Builds a row map of the saddle point system
Returns
true if successful, false otherwise
See Also
SetProblemMap , AddInterface , Mortar_Integrate
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

See Also
GetSPDRHSMatrix()

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:

  • Chooses a slave and a mortar side for each interface of not already chosen by the user in the interface class
  • Build a C0-continuous field of normals on the slave side
  • Makes special boundary modifications to the Lagrange multiplier shape functions
  • Performs integration of the mortar integrals on individual interfaces
  • Chooses degrees of freedom for the Lagrange multipliers based on the row map of the original uncoupled problem.
  • Builds a row map of the saddle point system
  • Assembles interface contributions to coupling matrices D and M
Returns
true if successful, false otherwise
See Also
SetProblemMap , AddInterface , Integrate_Interfaces
int MOERTEL::Manager::Ninterfaces ( )
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

Epetra_Map* MOERTEL::Manager::ProblemMap ( ) const
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.

See Also
SetProblemMap
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.

void MOERTEL::Manager::SetDimension ( MOERTEL::Manager::DimensionType  type)
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

Parameters
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

Parameters
map: Matrix of the original (uncoupled) problem
See Also
MakeSaddleProblem()
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()

Parameters
map: Row map of the original (uncoupled) problem
See Also
Mortar_Integrate()
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:

Teuchos::ParameterList params;
params.set("System","SaddleSystem"); // use a saddle point system of equations in solve
// or
params.set("System","SPDSystem"); // use a spd system with condensed Lagrange multipliers in solve
// choose solver package
params.set("Solver","Amesos"); // use solver package Amesos
// or
params.set("Solver","ML/Aztec"); // use solver packages ML and AztecOO
// argument sublist passed to and used for Amesos
// see Amesos documentation, all configured Amesos solvers can be used
// This is an example of possible parameters
Teuchos::ParameterList& amesosparams = params.sublist("Amesos");
amesosparams.set("Solver","Amesos_Klu"); // name of Amesos solver to use
amesosparams.set("PrintTiming",false);
amesosparams.set("PrintStatus",false);
amesosparams.set("UseTranspose",true);
// argument sublist passed to and used for ML
// see ML documentation, all configured parameters recognized by the
// ML_Epetra::ml_MultiLevelPreconditioner class can be used
// Moertel sets the ML default parameters first which are then overridden by
// this list
// This is an example configuration:
Teuchos::ParameterList& mlparams = params.sublist("ML");
ML_Epetra::SetDefaults("SA",mlparams);
mlparams.set("output",6);
mlparams.set("print unused",-2);
mlparams.set("increasing or decreasing","increasing");
mlparams.set("PDE equations",3);
mlparams.set("max levels",20);
mlparams.set("aggregation: type","Uncoupled");
mlparams.set("coarse: max size",2500);
mlparams.set("coarse: type","Amesos-KLU");
mlparams.set("smoother: type","MLS");
mlparams.set("smoother: MLS polynomial order",3);
mlparams.set("smoother: sweeps",1);
mlparams.set("smoother: pre or post","both");
mlparams.set("null space: type","pre-computed");
mlparams.set("null space: add default vectors",false);
int dimnullspace = 6;
int dimnsp = dimnullspace*nummyrows;
double* nsp = new double[dimnsp];
application_compute_nullspace(nsp,dimnsp);
mlparams.set("null space: dimension",dimnullspace);
// the user is responsible for freeing nsp after solve
mlparams.set("null space: vectors",nsp);
// argument sublist passed to and used for AztecOO
// see AztecOO documentation, all configured parameters recognized by the
// AztecOO class can be used
// Moertel sets the Aztec default parameters first which are then overridden by
// this list
// This is an example configuration:
Teuchos::ParameterList& aztecparams = params.sublist("Aztec");
aztecparams.set("AZ_solver","AZ_cg");
// if "AZ_user_precond" is specified, ML will be used as a preconditioner
// Any other aztec-buildin preconditioner can be used as well
aztecparams.set("AZ_precond","AZ_user_precond");
aztecparams.set("AZ_max_iter",500);
aztecparams.set("AZ_output",100);
aztecparams.set("AZ_tol",1.0e-08);
aztecparams.set("AZ_scaling","AZ_none");
Parameters
params: Parameter list containing solver options
Warning
When using ML and/or Aztec, one has to use "SPDSystem" as system matrix as the amg-preconditioned iterative method will currently fail on the saddle point system.
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.

Parameters
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.

Parameters
sol(out) : The solution
rhs(in) : The right hand side vector

References MOERTEL::Solver::Comm(), and MOERTEL::Solver::OutLevel().


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