Extension: Matrix functions
#include "dg/matrix/matrix.h"
dg::mat Namespace Reference

Classes

struct  BESSELI0
 \( f(x) = I_0 (x)\) with \(I_0\) the zeroth order modified Bessel function More...
 
struct  ConvertsToFunctionalButcherTableau
 Convert identifiers to their corresponding dg::mat::FunctionalButcherTableau. More...
 
struct  DirectSqrtCauchy
 Shortcut for \(b \approx \sqrt{A} x \) solve directly via sqrt Cauchy combined with PCG inversions. More...
 
struct  ExponentialERKStep
 Exponential Runge-Kutta fixed-step time-integration for \( \dot y = A y + g(t,y)\). More...
 
struct  ExponentialStep
 Exponential one step time-integration for \( \dot y = A y \). More...
 
struct  FunctionalButcherTableau
 Manage coefficients of a functional (extended) Butcher tableau. More...
 
struct  GAMMA0
 \( f(x) = \Gamma_0 (x) := I_0 (x) exp(x) \) with \(I_0\) the zeroth order modified Bessel function More...
 
struct  InvSqrtODE
 Right hand side of the square root ODE. More...
 
struct  MatrixFunction
 Computation of \( \vec x = f(A)\vec b\) for self-adjoint positive definite \( A\). More...
 
struct  MatrixSqrt
 Fast computation of \( \vec x = A^{\pm 1/2}\vec b\) for self-adjoint positive definite \(A\). More...
 
class  MCG
 EXPERIMENTAL Tridiagonalize \(A\) and approximate \(f(A)b \approx R f(\tilde T) e_1\) via CG algorithm. A is self-adjoint in the weights \( W\). More...
 
class  MCGFuncEigen
 EXPERIMENTAL Shortcut for \(x \approx f(A) b \) solve via exploiting first a Krylov projection achieved by the M-CG method and a matrix function computation via Eigen-decomposition. More...
 
class  PolCharge
 Various arbitary wavelength polarization charge operators of delta-f (df) and full-f (ff) More...
 
class  PolChargeN
 EXPERIMENTAL polarization solver class for N. More...
 
struct  SqrtCauchyInt
 Cauchy integral \( \sqrt{A} b= A\frac{ 2 K' \sqrt{m}}{\pi N} \sum_{j=1}^{N} ( w_j^2 I + A)^{-1} c_j d_j b \) More...
 
struct  TensorElliptic
 Matrix class that represents the arbitrary polarization operator. More...
 
class  TridiagInvD
 Compute the inverse of a general tridiagonal matrix. More...
 
class  TridiagInvDF
 USE THIS ONE Compute the inverse of a general tridiagonal matrix. The algorithm does not rely on the determinant. More...
 
class  TridiagInvHMGTI
 Compute the inverse of a general tridiagonal matrix. More...
 
class  UniversalLanczos
 Tridiagonalize \(A\) and approximate \(f(A)b \approx |b|_W V f(T) e_1\) via Lanczos algorithm. A is self-adjoint in the weights \( W\). More...
 

Enumerations

enum  func_tableau_identifier { EXPLICIT_EULER_1_1 , MIDPOINT_2_2 , CLASSIC_4_4 , HOCHBRUCK_3_3_4 }
 Identifiers for Butcher Tableaus. More...
 

Functions

template<class T >
phi1 (T x)
 \( f(x) = (\exp(x) - 1)/x\) More...
 
template<class T >
phi2 (T x)
 \( f(x) = (\exp(x) - x - 1)/x^2\) More...
 
template<class T >
phi3 (T x)
 \( f(x) = (\exp(x) - x^2/2 -x- 1)/x^3\) More...
 
template<class T >
phi4 (T x)
 \( f(x) = (\exp(x) - x^3/6 - x^2/2 -x- 1)/x^4\) More...
 
template<class UnaryOp >
auto make_FuncEigen_Te1 (UnaryOp f)
 Create a functor that uses Eigenvalue decomposition to compute \( f(T)\vec e_1 = E f(\Lambda) E^T \vec e_1 \) for symmetric tridiagonal T. More...
 
template<class value_type >
auto make_SqrtCauchy_Te1 (int exp, std::array< value_type, 2 > EVs, unsigned stepsCauchy)
 Create a functor that computes \( \sqrt{T^{\pm 1}} \vec e_1\) using SqrtCauchyInt. More...
 
template<class value_type >
auto make_SqrtCauchyEigen_Te1 (int exp, std::array< value_type, 2 > EVs, unsigned stepsCauchy)
 Create a functor that computes \( \sqrt{T^{\pm 1}} \vec e_1\) using either Eigen or SqrtCauchy solve based on whichever is fastest for given size. More...
 
template<class value_type >
auto make_SqrtODE_Te1 (int exp, std::string tableau, value_type rtol, value_type atol, unsigned &number)
 Create a functor that computes \( \sqrt{T^{\pm 1}} \vec e_1\) using ODE solve. More...
 
static auto make_Linear_Te1 (int exp)
 Create a functor that computes \( T^{\pm 1} \vec e_1\) directly. More...
 
template<class value_type , class ExplicitRHS >
auto make_directODESolve (ExplicitRHS &&ode, std::string tableau, value_type epsTimerel, value_type epsTimeabs, unsigned &number, value_type t0=0., value_type t1=1.)
 Operator that integrates an ODE from 0 to 1 with an adaptive ERK class as timestepper More...
 
template<class Matrix , class Preconditioner , class Container >
InvSqrtODE< Container > make_inv_sqrtodeCG (Matrix &A, const Preconditioner &P, const Container &weights, dg::get_value_type< Container > epsCG)
 Right hand side of the square root ODE. More...
 
template<class Matrix , class Container >
InvSqrtODE< Container > make_inv_sqrtodeTri (const Matrix &TH, const Container &copyable)
 
template<class MatrixType >
auto make_expode (MatrixType &A)
 Right hand side of the exponential ODE. More...
 
template<class MatrixType >
auto make_besselI0ode (MatrixType &A)
 Right hand side of the (zeroth order) modified Bessel function ODE, rewritten as a system of coupled first order ODEs: More...
 
template<class value_type >
value_type 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 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 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 > 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 > compute_extreme_EV (const cusp::dia_matrix< int, value_type, cusp::host_memory > &T)
 Compute extreme Eigenvalues of a symmetric tridiangular matrix. More...
 

Function Documentation

◆ make_inv_sqrtodeTri()

template<class Matrix , class Container >
InvSqrtODE< Container > dg::mat::make_inv_sqrtodeTri ( const Matrix &  TH,
const Container &  copyable 
)

◆ phi1()

template<class T >
T dg::mat::phi1 ( x)

\( f(x) = (\exp(x) - 1)/x\)

Accurate evaluation close to and at 0

Template Parameters
Tvalue type

◆ phi2()

template<class T >
T dg::mat::phi2 ( x)

\( f(x) = (\exp(x) - x - 1)/x^2\)

Accurate evaluation close to and at 0

Template Parameters
Tvalue type

◆ phi3()

template<class T >
T dg::mat::phi3 ( x)

\( f(x) = (\exp(x) - x^2/2 -x- 1)/x^3\)

Accurate evaluation close to and at 0

Template Parameters
Tvalue type

◆ phi4()

template<class T >
T dg::mat::phi4 ( x)

\( f(x) = (\exp(x) - x^3/6 - x^2/2 -x- 1)/x^4\)

Accurate evaluation close to and at 0

Template Parameters
Tvalue type