Amesos  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Amesos_Btf Class Reference

Amesos_Btf: Factors and solves a matrix after converting it to block triangular form. More...

#include <Amesos_BTF.h>

Inheritance diagram for Amesos_Btf:
Inheritance graph
[legend]
Collaboration diagram for Amesos_Btf:
Collaboration graph
[legend]

Public Member Functions

 Amesos_Btf (const Epetra_LinearProblem &LinearProblem)
 Amesos_Btf Constructor. More...
 
 ~Amesos_Btf (void)
 Amesos_Btf Destructor. More...
 
int SymbolicFactorization ()
 Performs SymbolicFactorization on the matrix A. More...
 
int NumericFactorization ()
 Performs NumericFactorization on the matrix A. More...
 
int Solve ()
 Solves A X = B (or AT X = B) More...
 
const Epetra_LinearProblemGetProblem () const
 Get a pointer to the Problem.
 
bool MatrixShapeOK () const
 Returns true if BTF can handle this matrix shape. More...
 
int SetUseTranspose (bool UseTranspose)
 SetUseTranpose(true) causes Solve() To compute A^T X = B. More...
 
bool UseTranspose () const
 Returns the current UseTranspose setting.
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix.
 
int SetParameters (Teuchos::ParameterList &ParameterList)
 Updates internal variables. More...
 
int NumSymbolicFact () const
 Returns the number of symbolic factorizations performed by this object.
 
int NumNumericFact () const
 Returns the number of numeric factorizations performed by this object.
 
int NumSolve () const
 Returns the number of solves performed by this object.
 
void PrintTiming ()
 Print timing information.
 
void PrintStatus ()
 Print information about the factorization and solution phases.
 
- Public Member Functions inherited from Amesos_BaseSolver
virtual ~Amesos_BaseSolver ()
 Destructor.
 
virtual void PrintStatus () const =0
 Prints status information about the current solver.
 
virtual void PrintTiming () const =0
 Prints timing information about the current solver.
 
virtual void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 Redefined from Teuchos::ParameterListAcceptor (Does Not Work)
 
virtual Teuchos::RCP
< Teuchos::ParameterList
getNonconstParameterList ()
 This is an empty stub.
 
virtual Teuchos::RCP
< Teuchos::ParameterList
unsetParameterList ()
 This is an empty stub.
 
virtual void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. (Does Not Work)
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 

Detailed Description

Amesos_Btf: Factors and solves a matrix after converting it to block triangular form.

Amesos_Btf:

Constructor & Destructor Documentation

Amesos_Btf::Amesos_Btf ( const Epetra_LinearProblem LinearProblem)

Amesos_Btf Constructor.

Creates an Amesos_Btf instance, using an Epetra_LinearProblem, passing in an already-defined Epetra_LinearProblem object.

Note: The operator in LinearProblem must be an Epetra_RowMatrix.

Amesos_Btf::~Amesos_Btf ( void  )

Amesos_Btf Destructor.

Completely deletes an Amesos_Btf object.

Member Function Documentation

bool Amesos_Btf::MatrixShapeOK ( ) const
virtual

Returns true if BTF can handle this matrix shape.

Returns true if the matrix shape is one that BTF can handle. BTF only works with square matrices.

Implements Amesos_BaseSolver.

int Amesos_Btf::NumericFactorization ( )
virtual

Performs NumericFactorization on the matrix A.

 \return Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

int Amesos_Btf::SetParameters ( Teuchos::ParameterList ParameterList)
virtual

Updates internal variables.

  <br \>Preconditions:<ul>
  <li>None.</li>
  </ul>

  <br \>Postconditions:<ul> 
  <li>Internal variables controlling the factorization and solve will
  be updated and take effect on all subsequent calls to NumericFactorization() 
  and Solve().</li>
  <li>All parameters whose value are to differ from the default values must 

be included in ParameterList. Parameters not specified in ParameterList revert to their default values.

Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

int Amesos_Btf::SetUseTranspose ( bool  UseTranspose)
inlinevirtual

SetUseTranpose(true) causes Solve() To compute A^T X = B.

Implements Amesos_BaseSolver.

References UseTranspose().

int Amesos_Btf::Solve ( )
virtual

Solves A X = B (or AT X = B)

Foreach block i:
   For each block j
     Compute x_i -= A_{i,j} x_j 
   Call Solve(x_i,b_i) 
   Broadcast x_i
 \return Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

int Amesos_Btf::SymbolicFactorization ( )
virtual

Performs SymbolicFactorization on the matrix A.

  • Compute an ordering which reduces the matrix to block upper triangular form.
  • Determine a task partitioning (static load balancing)
  • Redistribute the data based on the owner computes rule
  • Instantiates an Amesos solver for each of the diagonal blocks
  • Calls SymbolicFactorization() on each of the diagonal blocks
\return Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.


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