Extension: Matrix functions
#include "dg/matrix/matrix.h"
Inversion of tridiagonal matrices

\(T^{-1}\) More...

Collaboration diagram for Inversion of tridiagonal matrices:

Classes

class  dg::mat::TridiagInvHMGTI< ContainerType, DiaMatrix, CooMatrix >
 Compute the inverse of a general tridiagonal matrix. More...
 
class  dg::mat::TridiagInvDF< ContainerType, DiaMatrix, CooMatrix >
 USE THIS ONE Compute the inverse of a general tridiagonal matrix. The algorithm does not rely on the determinant. More...
 
class  dg::mat::TridiagInvD< ContainerType, DiaMatrix, CooMatrix >
 Compute the inverse of a general tridiagonal matrix. More...
 

Functions

template<class value_type >
value_type dg::mat::compute_Tinv_m1 (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T, unsigned size)
 Computes the value of \( (T^{-1})_{m1} = \langle \vec e_m, T^{-1}\vec e_1\rangle\) via a Thomas algorithm. More...
 
template<class value_type >
void dg::mat::compute_Tinv_y (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T, thrust::host_vector< value_type > &x, const thrust::host_vector< value_type > &y, value_type a=1., value_type d=0.)
 Computes the value of \( x = ((aT+dI)^{-1})y \) via Thomas algorithm. More...
 
template<class value_type >
void dg::mat::invert (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T, cusp::coo_matrix< int, value_type, cusp::host_memory > &Tinv)
 Invert a tridiagonal matrix. More...
 
template<class value_type >
cusp::coo_matrix< int, value_type, cusp::host_memory > dg::mat::invert (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T)
 Invert a tridiagonal matrix. More...
 
template<class value_type >
std::array< value_type, 2 > dg::mat::compute_extreme_EV (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T)
 Compute extreme Eigenvalues of a symmetric tridiangular matrix. More...
 

Detailed Description

\(T^{-1}\)

Function Documentation

◆ compute_extreme_EV()

template<class value_type >
std::array< value_type, 2 > dg::mat::compute_extreme_EV ( const cusp::dia_matrix< int, value_type, cusp::host_memory > &  T)

Compute extreme Eigenvalues of a symmetric tridiangular matrix.

dg::mat::UniversalLanczos lanczos( A.weights(), 20);
auto T = lanczos.tridiag( A, A.weights(), A.weights());
// EV[0] is the minimum, EV[1] the maximum Eigenvalue
Tridiagonalize and approximate via Lanczos algorithm. A is self-adjoint in the weights .
Definition: lanczos.h:68
std::array< value_type, 2 > compute_extreme_EV(const cusp::dia_matrix< int, value_type, cusp::host_memory > &T)
Compute extreme Eigenvalues of a symmetric tridiangular matrix.
Definition: tridiaginv.h:631
Template Parameters
value_typereal type
Parameters
Tsymmetric tridiangular matrix
Returns
{EVmin, EVmax}

◆ compute_Tinv_m1()

template<class value_type >
value_type dg::mat::compute_Tinv_m1 ( const cusp::dia_matrix< int, value_type, cusp::host_memory > &  T,
unsigned  size 
)

Computes the value of \( (T^{-1})_{m1} = \langle \vec e_m, T^{-1}\vec e_1\rangle\) via a Thomas algorithm.

Note
This is extremely fast (timings can bearly be measured; for size = 100 <1e-6s)
Template Parameters
value_typereal type
Parameters
TThe tridiagonal matrix
sizethe parameter m
Returns
\( (T^{-1})_{1m}\)

◆ compute_Tinv_y()

template<class value_type >
void dg::mat::compute_Tinv_y ( const cusp::dia_matrix< int, value_type, cusp::host_memory > &  T,
thrust::host_vector< value_type > &  x,
const thrust::host_vector< value_type > &  y,
value_type  a = 1.,
value_type  d = 0. 
)

Computes the value of \( x = ((aT+dI)^{-1})y \) via Thomas algorithm.

Note
This is extremely fast (timings can bearly be measured; for size = 100 <1e-6s)
Template Parameters
value_typereal type
Parameters
TThe tridiagonal matrix
xcontains the solution (resized if necessary)
ythe right hand side
aoptional scaling factor of T
doptional addition to diagonal of T

◆ invert() [1/2]

template<class value_type >
cusp::coo_matrix< int, value_type, cusp::host_memory > dg::mat::invert ( const cusp::dia_matrix< int, value_type, cusp::host_memory > &  T)

Invert a tridiagonal matrix.

This is a convenience shortcut for

return dg::TridiagInvDF( size)(T);
Template Parameters
value_typereal type
Parameters
T
Returns
Tinv

◆ invert() [2/2]

template<class value_type >
void dg::mat::invert ( const cusp::dia_matrix< int, value_type, cusp::host_memory > &  T,
cusp::coo_matrix< int, value_type, cusp::host_memory > &  Tinv 
)

Invert a tridiagonal matrix.

This is a convenience shortcut for

return dg::TridiagInvDF( size)(T, Tinv);
Template Parameters
value_typereal type
Parameters
T
Tinv(gets resized if necessary)