SparseDirectUMFPACK Class Reference
[Linear solver classesPreconditioners]

Inheritance diagram for SparseDirectUMFPACK:

Inheritance graph
[legend]

List of all members.

Classes

class  AdditionalData
class  ExcUMFPACKError

Public Member Functions

 SparseDirectUMFPACK ()
 ~SparseDirectUMFPACK ()
void initialize (const SparsityPattern &sparsity_pattern)
template<class Matrix >
void factorize (const Matrix &matrix)
template<class Matrix >
void initialize (const Matrix &matrix, const AdditionalData additional_data=AdditionalData())
void vmult (Vector< double > &, const Vector< double > &) const
void Tvmult (Vector< double > &, const Vector< double > &) const
void vmult_add (Vector< double > &, const Vector< double > &) const
void Tvmult_add (Vector< double > &, const Vector< double > &) const
void solve (Vector< double > &rhs_and_solution) const
template<class Matrix >
void solve (const Matrix &matrix, Vector< double > &rhs_and_solution)

Private Member Functions

void clear ()
template<typename number >
void sort_arrays (const SparseMatrix< number > &)
template<typename number >
void sort_arrays (const BlockSparseMatrix< number > &)

Private Attributes

void * symbolic_decomposition
void * numeric_decomposition
std::vector< long intAp
std::vector< long intAi
std::vector< doubleAx
std::vector< doublecontrol


Detailed Description

This class provides an interface to the sparse direct solver UMFPACK (see this link). UMFPACK is a set of routines for solving non-symmetric sparse linear systems, Ax=b, using the Unsymmetric-pattern MultiFrontal method and direct sparse LU factorization. Matrices may have symmetric or unsymmetrix sparsity patterns, and may have unsymmetric entries. The use of this class is explained in the step-22 and step-29 tutorial programs.

This matrix class implements the usual interface of preconditioners, that is a function initialize(const SparseMatrix<double>&matrix,const AdditionalData) for initalizing and the whole set of vmult() functions common to all matrices. Implemented here are only vmult() and vmult_add(), which perform multiplication with the inverse matrix. Furthermore, this class provides an older interface, consisting of the functions factorize() and solve(). Both interfaces are interchangeable.

Note:
This class only exists if support for UMFPACK was enabled during configure and if the UMFPACK library was configured. The steps to do this are explained in the deal.II ReadMe file. If you do nothing at the time you configure deal.II, then this class will simply not work.

UMFPACK has its own license, independent of that of deal.II. If you want to use the UMFPACK you have to accept that license. It is linked to from the deal.II ReadMe file. UMFPACK is included courtesy of its author, Timothy A. Davis.

Instantiations

There are instantiations of this class for SparseMatrix<double>, SparseMatrix<float>, BlockSparseMatrix<double>, and BlockSparseMatrix<float>.

See also:
UMFPACK
Author:
Wolfgang Bangerth, 2004

Constructor & Destructor Documentation

SparseDirectUMFPACK::SparseDirectUMFPACK (  ) 

Constructor. See the documentation of this class for the meaning of the parameters to this function.

SparseDirectUMFPACK::~SparseDirectUMFPACK (  ) 

Destructor.


Member Function Documentation

void SparseDirectUMFPACK::initialize ( const SparsityPattern sparsity_pattern  ) 

This function does nothing. It is only here to provide an interface that is consistent with that of the HSL MA27 and MA47 solver classes.

template<class Matrix >
void SparseDirectUMFPACK::factorize ( const Matrix &  matrix  )  [inline]

Factorize the matrix. This function may be called multiple times for different matrices, after the object of this class has been initialized for a certain sparsity pattern. You may therefore save some computing time if you want to invert several matrices with the same sparsity pattern. However, note that the bulk of the computing time is actually spent in the factorization, so this functionality may not always be of large benefit.

In contrast to the other direct solver classes, the initialisation method does nothing. Therefore initialise is not automatically called by this method, when the initialization step has not been performed yet.

This function copies the contents of the matrix into its own storage; the matrix can therefore be deleted after this operation, even if subsequent solves are required.

template<class Matrix >
void SparseDirectUMFPACK::initialize ( const Matrix &  matrix,
const AdditionalData  additional_data = AdditionalData() 
) [inline]

Initialize memory and call SparseDirectUMFPACK::factorize.

void SparseDirectUMFPACK::vmult ( Vector< double > &  ,
const Vector< double > &   
) const

Preconditioner interface function. Usually, given the source vector, this method returns an approximated solution of Ax = b. As this class provides a wrapper to a direct solver, here it is actually the exact solution (exact within the range of numerical accuracy of course).

void SparseDirectUMFPACK::Tvmult ( Vector< double > &  ,
const Vector< double > &   
) const

Not implemented but necessary for compiling.

void SparseDirectUMFPACK::vmult_add ( Vector< double > &  ,
const Vector< double > &   
) const

Same as vmult(), but adding to the previous solution. Not implemented yet.

void SparseDirectUMFPACK::Tvmult_add ( Vector< double > &  ,
const Vector< double > &   
) const

Not implemented but necessary for compiling.

void SparseDirectUMFPACK::solve ( Vector< double > &  rhs_and_solution  )  const

Solve for a certain right hand side vector. This function may be called multiple times for different right hand side vectors after the matrix has been factorized. This yields a big saving in computing time, since the actual solution is fast, compared to the factorization of the matrix.

The solution will be returned in place of the right hand side vector.

If the factorization has not happened before, strange things will happen. Note that we can't actually call the factorize() function from here if it has not yet been called, since we have no access to the actual matrix.

template<class Matrix >
void SparseDirectUMFPACK::solve ( const Matrix &  matrix,
Vector< double > &  rhs_and_solution 
) [inline]

Call the two functions factorize and solve in that order, i.e. perform the whole solution process for the given right hand side vector.

The solution will be returned in place of the right hand side vector.

void SparseDirectUMFPACK::clear (  )  [private]

Free all memory that hasn't been freed yet.

template<typename number >
void SparseDirectUMFPACK::sort_arrays ( const SparseMatrix< number > &   )  [inline, private]

Make sure that the arrays Ai and Ap are sorted in each row. UMFPACK wants it this way. We need to have two versions of this function, one for the usual SparseMatrix and one for the BlockSparseMatrix classes

template<typename number >
void SparseDirectUMFPACK::sort_arrays ( const BlockSparseMatrix< number > &   )  [inline, private]


Member Data Documentation

The UMFPACK routines allocate objects in which they store information about symbolic and numeric values of the decomposition. The actual data type of these objects is opaque, and only passed around as void pointers.

std::vector<long int> SparseDirectUMFPACK::Ap [private]

The arrays in which we store the data for the solver.

std::vector<long int> SparseDirectUMFPACK::Ai [private]

std::vector<double> SparseDirectUMFPACK::Ax [private]

std::vector<double> SparseDirectUMFPACK::control [private]

Control and info arrays for the solver routines.


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

deal.II documentation generated on Sat Aug 15 16:52:20 2009 by doxygen 1.5.9