![]() |
Reference documentation for deal.II version 8.1.0
|
#include <sparse_direct.h>
Classes | |
class | AdditionalData |
Public Types | |
typedef types::global_dof_index | size_type |
Public Member Functions | |
SparseDirectUMFPACK () | |
~SparseDirectUMFPACK () | |
DeclException2 (ExcUMFPACKError, char *, int,<< "UMFPACK routine "<< arg1<< " returned error status "<< arg2<< ". See the file <bundled/umfpack/UMFPACK/Include/umfpack.h>"<< " for a description of 'status codes'.") | |
Setting up a sparse factorization | |
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()) |
Functions that represent the inverse of a matrix | |
void | vmult (Vector< double > &dst, const Vector< double > &src) const |
void | vmult (BlockVector< double > &dst, const BlockVector< double > &src) const |
void | Tvmult (Vector< double > &dst, const Vector< double > &src) const |
void | Tvmult (BlockVector< double > &dst, const BlockVector< double > &src) const |
void | vmult_add (Vector< double > &dst, const Vector< double > &src) const |
void | Tvmult_add (Vector< double > &dst, const Vector< double > &src) const |
Functions that solve linear systems | |
void | solve (Vector< double > &rhs_and_solution, bool transpose=false) const |
void | solve (BlockVector< double > &rhs_and_solution, bool transpose=false) const |
template<class Matrix > | |
void | solve (const Matrix &matrix, Vector< double > &rhs_and_solution, bool transpose=false) |
template<class Matrix > | |
void | solve (const Matrix &matrix, BlockVector< double > &rhs_and_solution, bool transpose=false) |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
void | clear () |
template<typename number > | |
void | sort_arrays (const SparseMatrixEZ< number > &) |
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 int > | Ap |
std::vector< long int > | Ai |
std::vector< double > | Ax |
std::vector< double > | control |
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.
There are instantiations of this class for SparseMatrix<double>, SparseMatrix<float>, SparseMatrixEZ<float>, SparseMatrixEZ<double>, BlockSparseMatrix<double>, and BlockSparseMatrix<float>.
Definition at line 83 of file sparse_direct.h.
Declare type for container size.
Definition at line 89 of file sparse_direct.h.
SparseDirectUMFPACK::SparseDirectUMFPACK | ( | ) |
Constructor. See the documentation of this class for the meaning of the parameters to this function.
SparseDirectUMFPACK::~SparseDirectUMFPACK | ( | ) |
Destructor.
void SparseDirectUMFPACK::initialize | ( | const SparsityPattern & | sparsity_pattern | ) |
This function does nothing. It is only here to provide a interface consistent with other sparse direct solvers.
void SparseDirectUMFPACK::factorize | ( | const Matrix & | matrix | ) |
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.
void SparseDirectUMFPACK::initialize | ( | const Matrix & | matrix, |
const AdditionalData | additional_data = AdditionalData() |
||
) |
Initialize memory and call SparseDirectUMFPACK::factorize.
Preconditioner interface function. Usually, given the source vector, this method returns an approximate 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).
In other words, this function actually multiplies with the exact inverse of the matrix, .
void SparseDirectUMFPACK::vmult | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src | ||
) | const |
Same as before, but for block vectors.
Same as before, but uses the transpose of the matrix, i.e. this function multiplies with .
void SparseDirectUMFPACK::Tvmult | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src | ||
) | const |
Same as before, but for block vectors
Same as vmult(), but adding to the previous solution. Not implemented yet but necessary for compiling certain other classes.
Same as before, but uses the transpose of the matrix, i.e. this function multiplies with .
void SparseDirectUMFPACK::solve | ( | Vector< double > & | rhs_and_solution, |
bool | transpose = false |
||
) | 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.
If transpose
is set to true this function solves for the transpose of the matrix, i.e. .
void SparseDirectUMFPACK::solve | ( | BlockVector< double > & | rhs_and_solution, |
bool | transpose = false |
||
) | const |
Same as before, but for block vectors.
void SparseDirectUMFPACK::solve | ( | const Matrix & | matrix, |
Vector< double > & | rhs_and_solution, | ||
bool | transpose = false |
||
) |
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::solve | ( | const Matrix & | matrix, |
BlockVector< double > & | rhs_and_solution, | ||
bool | transpose = false |
||
) |
Same as before, but for block vectors.
SparseDirectUMFPACK::DeclException2 | ( | ExcUMFPACKError | , |
char * | , | ||
int | , | ||
<< "UMFPACK routine "<< arg1<< " returned error status "<< arg2<< ". See the file <bundled/umfpack/UMFPACK/Include/umfpack.h>"<< " for a description of 'status codes'." | |||
) |
One of the UMFPack routines threw an error. The error code is included in the output and can be looked up in the UMFPack user manual. The name of the routine is included for reference.
|
private |
Free all memory that hasn't been freed yet.
|
private |
Make sure that the arrays Ai and Ap are sorted in each row. UMFPACK wants it this way. We need to have three versions of this function, one for the usual SparseMatrix, one for the SparseMatrixEZ, and one for the BlockSparseMatrix classes
|
private |
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.
Definition at line 284 of file sparse_direct.h.
|
private |
The arrays in which we store the data for the solver.
Definition at line 310 of file sparse_direct.h.
|
private |
Control and info arrays for the solver routines.
Definition at line 317 of file sparse_direct.h.