Reference documentation for deal.II version 8.5.1
Private Attributes | List of all members
TridiagonalMatrix< number > Class Template Reference

#include <deal.II/lac/tridiagonal_matrix.h>

Public Types

Constructors
typedef types::global_dof_index size_type
 

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
 
Matrixnorms
number l1_norm () const
 
number linfty_norm () const
 
number frobenius_norm () const
 
number relative_symmetry_norm2 () 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
 
std::size_t memory_consumption () const
 

Private Attributes

std::vector< number > diagonal
 
std::vector< number > left
 
std::vector< number > right
 
bool is_symmetric
 
LAPACKSupport::State state
 

Detailed Description

template<typename number>
class TridiagonalMatrix< number >

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.

Author
Guido Kanschat, 2005, 2006

Definition at line 35 of file pointer_matrix.h.

Member Typedef Documentation

◆ size_type

template<typename number>
typedef types::global_dof_index TridiagonalMatrix< number >::size_type

Declare type for container size.

Definition at line 57 of file tridiagonal_matrix.h.

Constructor & Destructor Documentation

◆ TridiagonalMatrix()

template<typename number >
TridiagonalMatrix< number >::TridiagonalMatrix ( size_type  n = 0,
bool  symmetric = false 
)

Constructor generating an empty matrix of dimension n.

Definition at line 26 of file tridiagonal_matrix.cc.

Member Function Documentation

◆ reinit()

template<typename number >
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.

◆ m()

template<typename number>
size_type TridiagonalMatrix< number >::m ( ) const

Number of rows of this matrix. Note that the matrix is an m x m matrix.

◆ n()

template<typename number>
size_type TridiagonalMatrix< number >::n ( ) const

Number of columns of this matrix. Note that the matrix is an n x n matrix.

◆ all_zero()

template<typename number >
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 54 of file tridiagonal_matrix.cc.

◆ operator()() [1/2]

template<typename number>
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.

◆ operator()() [2/2]

template<typename number>
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.

Note
In case of symmetric storage technique, the entries (i,j) and (j,i) are identified and both exist. This must be taken into account if adding up is used for matrix assembling in order not to obtain doubled entries.

◆ vmult()

template<typename number >
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 78 of file tridiagonal_matrix.cc.

◆ vmult_add()

template<typename number >
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 132 of file tridiagonal_matrix.cc.

◆ Tvmult()

template<typename number >
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 142 of file tridiagonal_matrix.cc.

◆ Tvmult_add()

template<typename number >
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 186 of file tridiagonal_matrix.cc.

◆ matrix_scalar_product()

template<typename number >
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 196 of file tridiagonal_matrix.cc.

◆ matrix_norm_square()

template<typename number >
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 223 of file tridiagonal_matrix.cc.

◆ l1_norm()

template<typename number>
number TridiagonalMatrix< number >::l1_norm ( ) const

Return the \(l_1\)-norm of the matrix, i.e. \(|M|_1=max_{all columns j}\sum_{all rows i} |M_ij|\), (max. sum of columns). This is the natural matrix norm that is compatible to the \(l_1\)-norm for vectors, i.e. \(|Mv|_1\leq |M|_1 |v|_1\). (cf. Rannacher Numerik0)

◆ linfty_norm()

template<typename number>
number TridiagonalMatrix< number >::linfty_norm ( ) const

Return the \(l_\infty\)-norm of the matrix, i.e. \(|M|_\infty=\max_{all rows i}\sum_{all columns j} |M_{ij}|\), (max. sum of rows). This is the natural matrix norm that is compatible to the \(l_\infty\)-norm of vectors, i.e. \(|Mv|_\infty \leq |M|_\infty |v|_\infty\).

◆ frobenius_norm()

template<typename number>
number TridiagonalMatrix< number >::frobenius_norm ( ) const

The Frobenius norm of the matrix. Return value is the root of the square sum of all matrix entries.

◆ relative_symmetry_norm2()

template<typename number>
number TridiagonalMatrix< number >::relative_symmetry_norm2 ( ) const

Compute the relative norm of the skew-symmetric part. The return value is the Frobenius norm of the skew-symmetric part of the matrix divided by that of the matrix.

Main purpose of this function is to check, if a matrix is symmetric within a certain accuracy, or not.

◆ compute_eigenvalues()

template<typename number>
void TridiagonalMatrix< number >::compute_eigenvalues ( )

Compute the eigenvalues of the symmetric tridiagonal matrix.

Note
This function requires configuration of deal.II with LAPACK support. Additionally, the matrix must use symmetric storage technique.

◆ eigenvalue()

template<typename number >
number TridiagonalMatrix< number >::eigenvalue ( const size_type  i) const

After calling compute_eigenvalues(), you can access each eigenvalue here.

Definition at line 252 of file tridiagonal_matrix.cc.

◆ print()

template<typename number>
template<class OutputStream >
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.

◆ memory_consumption()

template<typename number>
std::size_t TridiagonalMatrix< number >::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Member Data Documentation

◆ diagonal

template<typename number>
std::vector<number> TridiagonalMatrix< number >::diagonal
private

The diagonal entries.

Definition at line 260 of file tridiagonal_matrix.h.

◆ left

template<typename number>
std::vector<number> TridiagonalMatrix< number >::left
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 270 of file tridiagonal_matrix.h.

◆ right

template<typename number>
std::vector<number> TridiagonalMatrix< number >::right
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 276 of file tridiagonal_matrix.h.

◆ is_symmetric

template<typename number>
bool TridiagonalMatrix< number >::is_symmetric
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 282 of file tridiagonal_matrix.h.

◆ state

template<typename number>
LAPACKSupport::State TridiagonalMatrix< number >::state
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 290 of file tridiagonal_matrix.h.


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