Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/tridiagonal_matrix.h>
Public Types | |
Constructors | |
using | size_type = types::global_dof_index |
Public Member Functions | |
Constructors and initialization | |
TridiagonalMatrix (size_type n=0, bool symmetric=false) | |
void | reinit (size_type n, bool symmetric=false) |
Non-modifying operators | |
size_type | m () const |
size_type | n () const |
bool | all_zero () const |
Element access | |
number | operator() (size_type i, size_type j) const |
number & | operator() (size_type i, size_type j) |
Multiplications with vectors | |
void | vmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
void | vmult_add (Vector< number > &w, const Vector< number > &v) const |
void | Tvmult (Vector< number > &w, const Vector< number > &v, const bool adding=false) const |
void | Tvmult_add (Vector< number > &w, const Vector< number > &v) const |
number | matrix_scalar_product (const Vector< number > &u, const Vector< number > &v) const |
number | matrix_norm_square (const Vector< number > &v) const |
LAPACK operations | |
void | compute_eigenvalues () |
number | eigenvalue (const size_type i) const |
Miscellanea | |
template<class OutputStream > | |
void | print (OutputStream &s, const unsigned int width=5, const unsigned int precision=2) const |
Private Attributes | |
std::vector< number > | diagonal |
std::vector< number > | left |
std::vector< number > | right |
bool | is_symmetric |
LAPACKSupport::State | state |
A quadratic tridiagonal matrix. That is, a matrix where all entries are zero, except the diagonal and the entries left and right of it.
The matrix has an additional symmetric mode, in which case only the upper triangle of the matrix is stored and mirrored to the lower one for matrix vector operations.
Definition at line 44 of file pointer_matrix.h.
using TridiagonalMatrix< number >::size_type = types::global_dof_index |
Declare type for container size.
Definition at line 60 of file tridiagonal_matrix.h.
TridiagonalMatrix< number >::TridiagonalMatrix | ( | size_type | n = 0 , |
bool | symmetric = false |
||
) |
Constructor generating an empty matrix of dimension n
.
Definition at line 28 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::reinit | ( | size_type | n, |
bool | symmetric = false |
||
) |
Reinitialize the matrix to a new size and reset all entries to zero. The symmetry properties may be set as well.
Definition at line 40 of file tridiagonal_matrix.cc.
size_type TridiagonalMatrix< number >::m | ( | ) | const |
Number of rows of this matrix. Note that the matrix is an m x m matrix.
size_type TridiagonalMatrix< number >::n | ( | ) | const |
Number of columns of this matrix. Note that the matrix is an n x n matrix.
bool TridiagonalMatrix< number >::all_zero | ( | ) | const |
Return whether the matrix contains only elements with value zero. This function is mainly for internal consistency checks and should seldom be used when not in debug mode since it uses quite some time.
Definition at line 53 of file tridiagonal_matrix.cc.
number TridiagonalMatrix< number >::operator() | ( | size_type | i, |
size_type | j | ||
) | const |
Read-only access to a value. This is restricted to the case where |i-j| <= 1.
number& TridiagonalMatrix< number >::operator() | ( | size_type | i, |
size_type | j | ||
) |
Read-write access to a value. This is restricted to the case where |i-j| <= 1.
void TridiagonalMatrix< number >::vmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false |
||
) | const |
Matrix-vector-multiplication. Multiplies v
from the right and stores the result in w
.
If the optional parameter adding
is true
, the result is added to w
.
Source and destination must not be the same vector.
Definition at line 81 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::vmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v | ||
) | const |
Adding Matrix-vector-multiplication. Same as vmult() with parameter adding=true
, but widely used in deal.II
classes.
Source and destination must not be the same vector.
Definition at line 135 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::Tvmult | ( | Vector< number > & | w, |
const Vector< number > & | v, | ||
const bool | adding = false |
||
) | const |
Transpose matrix-vector-multiplication. Multiplies vT
from the left and stores the result in w
.
If the optional parameter adding
is true
, the result is added to w
.
Source and destination must not be the same vector.
Definition at line 145 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::Tvmult_add | ( | Vector< number > & | w, |
const Vector< number > & | v | ||
) | const |
Adding transpose matrix-vector-multiplication. Same as Tvmult() with parameter adding=true
, but widely used in deal.II
classes.
Source and destination must not be the same vector.
Definition at line 190 of file tridiagonal_matrix.cc.
number TridiagonalMatrix< number >::matrix_scalar_product | ( | const Vector< number > & | u, |
const Vector< number > & | v | ||
) | const |
Build the matrix scalar product u^T M v
. This function is mostly useful when building the cellwise scalar product of two functions in the finite element context.
Definition at line 200 of file tridiagonal_matrix.cc.
number TridiagonalMatrix< number >::matrix_norm_square | ( | const Vector< number > & | v | ) | const |
Return the square of the norm of the vector v
with respect to the norm induced by this matrix, i.e. (v,Mv). This is useful, e.g. in the finite element context, where the L2 norm of a function equals the matrix norm with respect to the mass matrix of the vector representing the nodal values of the finite element function.
Obviously, the matrix needs to be quadratic for this operation.
Definition at line 227 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::compute_eigenvalues | ( | ) |
Compute the eigenvalues of the symmetric tridiagonal matrix.
Definition at line 236 of file tridiagonal_matrix.cc.
number TridiagonalMatrix< number >::eigenvalue | ( | const size_type | i | ) | const |
After calling compute_eigenvalues(), you can access each eigenvalue here.
Definition at line 264 of file tridiagonal_matrix.cc.
void TridiagonalMatrix< number >::print | ( | OutputStream & | s, |
const unsigned int | width = 5 , |
||
const unsigned int | precision = 2 |
||
) | const |
Output of the matrix in user-defined format.
|
private |
The diagonal entries.
Definition at line 235 of file tridiagonal_matrix.h.
|
private |
The entries left of the diagonal. The entry with index zero is always zero, since the first row has no entry left of the diagonal. Therefore, the length of this vector is the same as that of diagonal.
The length of this vector is zero for symmetric storage. In this case, the second element of left is identified with the first element of right.
Definition at line 245 of file tridiagonal_matrix.h.
|
private |
The entries right of the diagonal. The last entry is always zero, since the last row has no entry right of the diagonal. Therefore, the length of this vector is the same as that of diagonal.
Definition at line 251 of file tridiagonal_matrix.h.
|
private |
If this flag is true, only the entries to the right of the diagonal are stored and the matrix is assumed symmetric.
Definition at line 257 of file tridiagonal_matrix.h.
|
private |
The state of the matrix. Normally, the state is LAPACKSupport::matrix, indicating that the object can be used for regular matrix operations.
See explanation of this data type for details.
Definition at line 265 of file tridiagonal_matrix.h.