Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/lapack_full_matrix.h>
Public Types | |
typedef types::global_dof_index | size_type |
Public Types inherited from TransposeTable< number > | |
typedef TableBase< 2, number >::size_type | size_type |
Public Types inherited from TableBase< 2, number > | |
typedef AlignedVector< number >::size_type | size_type |
Public Member Functions | |
LAPACKFullMatrix (const size_type size=0) | |
LAPACKFullMatrix (const size_type rows, const size_type cols) | |
LAPACKFullMatrix (const LAPACKFullMatrix &) | |
LAPACKFullMatrix< number > & | operator= (const LAPACKFullMatrix< number > &) |
template<typename number2 > | |
LAPACKFullMatrix< number > & | operator= (const FullMatrix< number2 > &) |
template<typename number2 > | |
LAPACKFullMatrix< number > & | operator= (const SparseMatrix< number2 > &) |
LAPACKFullMatrix< number > & | operator= (const double d) |
LAPACKFullMatrix< number > & | operator*= (const number factor) |
LAPACKFullMatrix< number > & | operator/= (const number factor) |
template<typename MatrixType > | |
void | copy_from (const MatrixType &) |
void | reinit (const size_type size) |
void | reinit (const size_type rows, const size_type cols) |
unsigned int | m () const |
unsigned int | n () const |
template<typename MatrixType > | |
void | fill (const MatrixType &src, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0, const number factor=1., const bool transpose=false) |
template<typename number2 > | |
void | vmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
void | vmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
template<typename number2 > | |
void | vmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
void | vmult_add (Vector< number > &w, const Vector< number > &v) const |
template<typename number2 > | |
void | Tvmult (Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const |
void | Tvmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
template<typename number2 > | |
void | Tvmult_add (Vector< number2 > &w, const Vector< number2 > &v) const |
void | Tvmult_add (Vector< number > &w, const Vector< number > &v) const |
void | mmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | mmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | Tmmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | Tmmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | mTmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | mTmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | TmTmult (LAPACKFullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | TmTmult (FullMatrix< number > &C, const LAPACKFullMatrix< number > &B, const bool adding=false) const |
void | compute_lu_factorization () |
void | invert () |
void | apply_lu_factorization (Vector< number > &v, const bool transposed) const |
void | apply_lu_factorization (LAPACKFullMatrix< number > &B, const bool transposed) const |
void | compute_eigenvalues (const bool right_eigenvectors=false, const bool left_eigenvectors=false) |
void | compute_eigenvalues_symmetric (const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, FullMatrix< number > &eigenvectors) |
void | compute_generalized_eigenvalues_symmetric (LAPACKFullMatrix< number > &B, const number lower_bound, const number upper_bound, const number abs_accuracy, Vector< number > &eigenvalues, std::vector< Vector< number > > &eigenvectors, const int itype=1) |
void | compute_generalized_eigenvalues_symmetric (LAPACKFullMatrix< number > &B, std::vector< Vector< number > > &eigenvectors, const int itype=1) |
void | compute_svd () |
void | compute_inverse_svd (const double threshold=0.) |
std::complex< number > | eigenvalue (const size_type i) const |
number | singular_value (const size_type i) const |
void | print_formatted (std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const double threshold=0.) const |
Public Member Functions inherited from TransposeTable< number > | |
TransposeTable () | |
TransposeTable (const unsigned int size1, const unsigned int size2) | |
void | reinit (const unsigned int size1, const unsigned int size2, const bool omit_default_initialization=false) |
AlignedVector< number >::const_reference | operator() (const unsigned int i, const unsigned int j) const |
AlignedVector< number >::reference | operator() (const unsigned int i, const unsigned int j) |
unsigned int | n_rows () const |
unsigned int | n_cols () const |
Public Member Functions inherited from TableBase< 2, number > | |
TableBase () | |
TableBase (const TableIndices< N > &sizes) | |
TableBase (const TableIndices< N > &sizes, InputIterator entries, const bool C_style_indexing=true) | |
TableBase (const TableBase< N, number > &src) | |
TableBase (const TableBase< N, T2 > &src) | |
TableBase (TableBase< N, number > &&src) | |
~TableBase () | |
TableBase< N, number > & | operator= (const TableBase< N, number > &src) |
TableBase< N, number > & | operator= (const TableBase< N, T2 > &src) |
TableBase< N, number > & | operator= (TableBase< N, number > &&src) |
bool | operator== (const TableBase< N, number > &T2) const |
void | reset_values () |
void | reinit (const TableIndices< N > &new_size, const bool omit_default_initialization=false) |
unsigned int | size (const unsigned int i) const |
const TableIndices< N > & | size () const |
size_type | n_elements () const |
bool | empty () const |
void | fill (InputIterator entries, const bool C_style_indexing=true) |
void | fill (const number &value) |
AlignedVector< number >::reference | operator() (const TableIndices< N > &indices) |
AlignedVector< number >::const_reference | operator() (const TableIndices< N > &indices) const |
void | swap (TableBase< N, number > &v) |
std::size_t | memory_consumption () const |
void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Attributes | |
LAPACKSupport::State | state |
LAPACKSupport::Properties | properties |
std::vector< number > | work |
std::vector< int > | ipiv |
std::vector< number > | inv_work |
std::vector< number > | wr |
std::vector< number > | wi |
std::vector< number > | vl |
std::vector< number > | vr |
std_cxx11::shared_ptr< LAPACKFullMatrix< number > > | svd_u |
std_cxx11::shared_ptr< LAPACKFullMatrix< number > > | svd_vt |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Protected Member Functions inherited from TransposeTable< number > | |
AlignedVector< number >::reference | el (const unsigned int i, const unsigned int j) |
AlignedVector< number >::const_reference | el (const unsigned int i, const unsigned int j) const |
Protected Member Functions inherited from TableBase< 2, number > | |
size_type | position (const TableIndices< N > &indices) const |
AlignedVector< number >::reference | el (const TableIndices< N > &indices) |
AlignedVector< number >::const_reference | el (const TableIndices< N > &indices) const |
Protected Attributes inherited from TableBase< 2, number > | |
AlignedVector< number > | values |
TableIndices< N > | table_size |
A variant of FullMatrix using LAPACK functions wherever possible. In order to do this, the matrix is stored in transposed order. The element access functions hide this fact by reverting the transposition.
Definition at line 36 of file full_matrix.h.
typedef types::global_dof_index LAPACKFullMatrix< number >::size_type |
Declare type for container size.
Definition at line 59 of file lapack_full_matrix.h.
|
explicit |
Constructor. Initialize the matrix as a square matrix with dimension n
.
In order to avoid the implicit conversion of integers and other types to a matrix, this constructor is declared explicit
.
By default, no memory is allocated.
Definition at line 33 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number >::LAPACKFullMatrix | ( | const size_type | rows, |
const size_type | cols | ||
) |
Constructor. Initialize the matrix as a rectangular matrix.
Definition at line 41 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number >::LAPACKFullMatrix | ( | const LAPACKFullMatrix< number > & | M | ) |
Copy constructor. This constructor does a deep copy of the matrix. Therefore, it poses a possible efficiency problem, if for example, function arguments are passed by value rather than by reference. Unfortunately, we can't mark this copy constructor explicit
, since that prevents the use of this class in containers, such as std::vector
. The responsibility to check performance of programs must therefore remain with the user of this class.
Definition at line 50 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const LAPACKFullMatrix< number > & | M | ) |
Assignment operator.
Definition at line 59 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const FullMatrix< number2 > & | M | ) |
Assignment operator from a regular FullMatrix.
Definition at line 89 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const SparseMatrix< number2 > & | M | ) |
Assignment operator from a regular SparseMatrix.
Definition at line 105 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator= | ( | const double | d | ) |
This operator assigns a scalar to a matrix. To avoid confusion with constructors, zero is the only value allowed for d
Definition at line 120 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator*= | ( | const number | factor | ) |
This operator multiplies all entries by a fixed factor.
Definition at line 135 of file lapack_full_matrix.cc.
LAPACKFullMatrix< number > & LAPACKFullMatrix< number >::operator/= | ( | const number | factor | ) |
This operator divides all entries by a fixed factor.
Definition at line 151 of file lapack_full_matrix.cc.
|
inline |
Assignment from different matrix classes, performing the usual conversion to the transposed format expected by LAPACK. This assignment operator uses iterators of the typename MatrixType. Therefore, sparse matrices are possible sources.
Definition at line 723 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::reinit | ( | const size_type | size | ) |
Regenerate the current matrix by one that has the same properties as if it were created by the constructor of this class with the same argument list as this present function.
Definition at line 69 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::reinit | ( | const size_type | rows, |
const size_type | cols | ||
) |
Regenerate the current matrix by one that has the same properties as if it were created by the constructor of this class with the same argument list as this present function.
Definition at line 78 of file lapack_full_matrix.cc.
|
inline |
Return the dimension of the codomain (or range) space.
Definition at line 707 of file lapack_full_matrix.h.
|
inline |
Return the dimension of the domain space.
Definition at line 715 of file lapack_full_matrix.h.
|
inline |
Fill rectangular block.
A rectangular block of the matrix src
is copied into this
. The upper left corner of the block being copied is (src_offset_i,src_offset_j)
. The upper left corner of the copied block is (dst_offset_i,dst_offset_j)
. The size of the rectangular block being copied is the maximum size possible, determined either by the size of this
or src
.
The final two arguments allow to enter a multiple of the source or its transpose.
Definition at line 746 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v, | ||
const bool | adding = false |
||
) | const |
Matrix-vector-multiplication.
Depending on previous transformations recorded in state, the result of this function is one of
The optional parameter adding
determines, whether the result is stored in w
or added to w
.
if (adding) w += A*v
if (!adding) w = A*v
number2
only exists for compile-time compatibility with FullMatrix. Only the case number2
= number
is implemented due to limitations in the underlying LAPACK interface. All other variants throw an error upon invocation. Definition at line 779 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false |
||
) | const |
Specialization of above function for compatible Vector::value_type.
Definition at line 170 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::vmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v | ||
) | const |
Adding Matrix-vector-multiplication. w += A*v
See the documentation of vmult() for details on the implementation.
Definition at line 792 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::vmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v | ||
) | const |
Specialization of above function for compatible Vector::value_type.
Definition at line 288 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tvmult | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v, | ||
const bool | adding = false |
||
) | const |
Transpose matrix-vector-multiplication.
The optional parameter adding
determines, whether the result is stored in w
or added to w
.
if (adding) w += AT*v
if (!adding) w = AT*v
See the documentation of vmult() for details on the implementation.
Definition at line 804 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::Tvmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false |
||
) | const |
Specialization of above function for compatible Vector::value_type.
Definition at line 228 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tvmult_add | ( | Vector< number2 > & | w, |
const Vector< number2 > & | v | ||
) | const |
Adding transpose matrix-vector-multiplication. w += AT*v
See the documentation of vmult() for details on the implementation.
Definition at line 817 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::Tvmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v | ||
) | const |
Specialization of above function for compatible Vector::value_type.
Definition at line 297 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += A*B
if (!adding) C = A*B
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm.
Definition at line 306 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 329 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tmmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of this
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += AT*B
if (!adding) C = AT*B
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm.
Definition at line 354 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::Tmmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 377 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mTmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of B
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += A*BT
if (!adding) C = A*BT
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm.
Definition at line 401 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::mTmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 425 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::TmTmult | ( | LAPACKFullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Matrix-matrix-multiplication using transpose of this
and B
.
The optional parameter adding
determines, whether the result is stored in C
or added to C
.
if (adding) C += AT*BT
if (!adding) C = AT*BT
Assumes that A
and B
have compatible sizes and that C
already has the right size.
This function uses the BLAS function Xgemm.
Definition at line 449 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::TmTmult | ( | FullMatrix< number > & | C, |
const LAPACKFullMatrix< number > & | B, | ||
const bool | adding = false |
||
) | const |
Same as before, but stores the result in a FullMatrix, not in a LAPACKFullMatrix.
Definition at line 472 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_lu_factorization | ( | ) |
Compute the LU factorization of the matrix using LAPACK function Xgetrf.
Definition at line 496 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::invert | ( | ) |
Invert the matrix by first computing an LU factorization with the LAPACK function Xgetrf and then building the actual inverse using Xgetri.
Definition at line 585 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::apply_lu_factorization | ( | Vector< number > & | v, |
const bool | transposed | ||
) | const |
Solve the linear system with right hand side given by applying forward/backward substitution to the previously computed LU factorization. Uses LAPACK function Xgetrs.
The flag transposed indicates whether the solution of the transposed system is to be performed.
Definition at line 617 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::apply_lu_factorization | ( | LAPACKFullMatrix< number > & | B, |
const bool | transposed | ||
) | const |
Solve the linear system with multiple right hand sides (as many as there are columns in the matrix b) given by applying forward/backward substitution to the previously computed LU factorization. Uses LAPACK function Xgetrs.
The flag transposed indicates whether the solution of the transposed system is to be performed.
Definition at line 639 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_eigenvalues | ( | const bool | right_eigenvectors = false , |
const bool | left_eigenvectors = false |
||
) |
Compute eigenvalues of the matrix. After this routine has been called, eigenvalues can be retrieved using the eigenvalue() function. The matrix itself will be LAPACKSupport::unusable after this operation.
The optional arguments allow to compute left and right eigenvectors as well.
Note that the function does not return the computed eigenvalues right away since that involves copying data around between the output arrays of the LAPACK functions and any return array. This is often unnecessary since one may not be interested in all eigenvalues at once, but for example only the extreme ones. In that case, it is cheaper to just have this function compute the eigenvalues and have a separate function that returns whatever eigenvalue is requested.
Definition at line 661 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_eigenvalues_symmetric | ( | const number | lower_bound, |
const number | upper_bound, | ||
const number | abs_accuracy, | ||
Vector< number > & | eigenvalues, | ||
FullMatrix< number > & | eigenvectors | ||
) |
Compute eigenvalues and eigenvectors of a real symmetric matrix. Only eigenvalues in the interval (lower_bound, upper_bound] are computed with the absolute tolerance abs_accuracy. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a,b] of width less than or equal to abs_accuracy + eps * max( |a|,|b| ), where eps is the machine precision. If abs_accuracy is less than or equal to zero, then eps*|t| will be used in its place, where |t| is the 1-norm of the tridiagonal matrix obtained by reducing A to tridiagonal form. Eigenvalues will be computed most accurately when abs_accuracy is set to twice the underflow threshold, not zero. After this routine has been called, all eigenvalues in (lower_bound, upper_bound] will be stored in eigenvalues and the corresponding eigenvectors will be stored in the columns of eigenvectors, whose dimension is set accordingly.
Definition at line 722 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_generalized_eigenvalues_symmetric | ( | LAPACKFullMatrix< number > & | B, |
const number | lower_bound, | ||
const number | upper_bound, | ||
const number | abs_accuracy, | ||
Vector< number > & | eigenvalues, | ||
std::vector< Vector< number > > & | eigenvectors, | ||
const int | itype = 1 |
||
) |
Compute generalized eigenvalues and eigenvectors of a real generalized symmetric eigenproblem of the form itype = 1: \(Ax=\lambda B x\) itype = 2: \(ABx=\lambda x\) itype = 3: \(BAx=\lambda x\), where A is this matrix. A and B are assumed to be symmetric, and B has to be positive definite. Only eigenvalues in the interval (lower_bound, upper_bound] are computed with the absolute tolerance abs_accuracy. An approximate eigenvalue is accepted as converged when it is determined to lie in an interval [a,b] of width less than or equal to abs_accuracy + eps * max( |a|,|b| ), where eps is the machine precision. If abs_accuracy is less than or equal to zero, then eps*|t| will be used in its place, where |t| is the 1-norm of the tridiagonal matrix obtained by reducing A to tridiagonal form. Eigenvalues will be computed most accurately when abs_accuracy is set to twice the underflow threshold, not zero. After this routine has been called, all eigenvalues in (lower_bound, upper_bound] will be stored in eigenvalues and the corresponding eigenvectors will be stored in eigenvectors, whose dimension is set accordingly.
Definition at line 808 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_generalized_eigenvalues_symmetric | ( | LAPACKFullMatrix< number > & | B, |
std::vector< Vector< number > > & | eigenvectors, | ||
const int | itype = 1 |
||
) |
Same as the other compute_generalized_eigenvalues_symmetric function except that all eigenvalues are computed and the tolerance is set automatically. Note that this function does not return the computed eigenvalues right away since that involves copying data around between the output arrays of the LAPACK functions and any return array. This is often unnecessary since one may not be interested in all eigenvalues at once, but for example only the extreme ones. In that case, it is cheaper to just have this function compute the eigenvalues and have a separate function that returns whatever eigenvalue is requested. Eigenvalues can be retrieved using the eigenvalue() function. The number of computed eigenvectors is equal to eigenvectors.size()
Definition at line 900 of file lapack_full_matrix.cc.
void LAPACKFullMatrix< number >::compute_svd | ( | ) |
void LAPACKFullMatrix< number >::compute_inverse_svd | ( | const double | threshold = 0. | ) |
Compute the inverse of the matrix by singular value decomposition.
Requires that state is either LAPACKSupport::matrix or LAPACKSupport::svd. In the first case, this function calls compute_svd(). After this function, the object will have the state LAPACKSupport::inverse_svd.
For a singular value decomposition, the inverse is simply computed by replacing all singular values by their reciprocal values. If the matrix does not have maximal rank, singular values 0 are not touched, thus computing the minimal norm right inverse of the matrix.
The parameter threshold
determines, when a singular value should be considered zero. It is the ratio of the smallest to the largest nonzero singular value smax. Thus, the inverses of all singular values less than smax/threshold
will be set to zero.
Definition at line 564 of file lapack_full_matrix.cc.
|
inline |
Retrieve eigenvalue after compute_eigenvalues() was called.
Definition at line 828 of file lapack_full_matrix.h.
|
inline |
Retrieve singular values after compute_svd() or compute_inverse_svd() was called.
Definition at line 841 of file lapack_full_matrix.h.
void LAPACKFullMatrix< number >::print_formatted | ( | std::ostream & | out, |
const unsigned int | precision = 3 , |
||
const bool | scientific = true , |
||
const unsigned int | width = 0 , |
||
const char * | zero_string = " " , |
||
const double | denominator = 1. , |
||
const double | threshold = 0. |
||
) | const |
Print the matrix and allow formatting of entries.
The parameters allow for a flexible setting of the output format:
precision
denotes the number of trailing digits.scientific
is used to determine the number format, where scientific
= false
means fixed point notation.width
denotes the with of each column. A zero entry for width
makes the function compute a width, but it may be changed to a positive value, if output is crude.zero_string
specifies a string printed for zero entries.denominator
Multiply the whole matrix by this common denominator to get nicer numbers.threshold
: all entries with absolute value smaller than this are considered zero.Definition at line 975 of file lapack_full_matrix.cc.
|
private |
Since LAPACK operations notoriously change the meaning of the matrix entries, we record the current state after the last operation here.
Definition at line 614 of file lapack_full_matrix.h.
|
private |
Additional properties of the matrix which may help to select more efficient LAPACK functions.
Definition at line 620 of file lapack_full_matrix.h.
|
mutableprivate |
The working array used for some LAPACK functions.
Definition at line 625 of file lapack_full_matrix.h.
|
private |
The vector storing the permutations applied for pivoting in the LU- factorization.
Also used as the scratch array IWORK for LAPACK functions needing it.
Definition at line 633 of file lapack_full_matrix.h.
|
private |
Workspace for calculating the inverse matrix from an LU factorization.
Definition at line 638 of file lapack_full_matrix.h.
|
private |
Real parts of eigenvalues or the singular values. Filled by compute_eigenvalues() or compute_svd().
Definition at line 644 of file lapack_full_matrix.h.
|
private |
Imaginary parts of eigenvalues. Filled by compute_eigenvalues.
Definition at line 649 of file lapack_full_matrix.h.
|
private |
Space where left eigenvectors can be stored.
Definition at line 654 of file lapack_full_matrix.h.
|
private |
Space where right eigenvectors can be stored.
Definition at line 659 of file lapack_full_matrix.h.
|
private |
The matrix U in the singular value decomposition USVT.
Definition at line 665 of file lapack_full_matrix.h.
|
private |
The matrix VT in the singular value decomposition USVT.
Definition at line 671 of file lapack_full_matrix.h.