\( x \approx f(T)b \)
More...
|
template<class UnaryOp > |
auto | dg::mat::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 | dg::mat::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 | dg::mat::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 | dg::mat::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 | dg::mat::make_Linear_Te1 (int exp) |
| Create a functor that computes \( T^{\pm 1} \vec e_1\) directly. More...
|
|
\( x \approx f(T)b \)
approximation
◆ make_FuncEigen_Te1()
template<class UnaryOp >
auto dg::mat::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
.
- Note
- This is a general purpose solution. Very fast for small sizes (<40) of T, but scales badly for larger sizes. Use more specialized solutions if the number of iterations becomes high
- Parameters
-
f | the matrix function (e.g. dg::SQRT<double> or dg::EXP<double>) |
- Returns
- an operator to use in
UniversalLanczos
solve method
- See also
UniversalLanczos
◆ make_Linear_Te1()
static auto dg::mat::make_Linear_Te1 |
( |
int |
exp | ) |
|
|
inlinestatic |
Create a functor that computes \( T^{\pm 1} \vec e_1\) directly.
- Parameters
-
exp | exponent if +1 compute \( Te_1\), if -1 compute \( T^{-1} e_1\) |
- Returns
- an operator to use in
UniversalLanczos
solve method
- See also
UniversalLanczos
◆ make_SqrtCauchy_Te1()
template<class value_type >
auto dg::mat::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.
- Note
- The Eigenvalues can be estimated from a few lanczos iterations (which is at least more reliable than doing it semi-analytically)
auto T = lanczos.tridiag( A, A.weights(), A.weights());
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
auto make_SqrtCauchy_Te1(int exp, std::array< value_type, 2 > EVs, unsigned stepsCauchy)
Create a functor that computes using SqrtCauchyInt.
Definition: matrixfunction.h:92
- Parameters
-
exp | exponent if +1 compute \( \sqrt(T)\), if -1 compute \( 1/\sqrt(T)\) |
EVs | {minimum Eigenvalue of A, maximum Eigenvalue of A} |
stepsCauchy | iterations of cauchy integral |
- Returns
- an operator to use in
UniversalLanczos
solve method
- See also
SqrtCauchyInt
UniversalLanczos
◆ make_SqrtCauchyEigen_Te1()
template<class value_type >
auto dg::mat::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.
- Note
- The Eigenvalues can be estimated from a few lanczos iterations (which is at least more reliable than doing it semi-analytically)
auto T = lanczos.tridiag( A, A.weights(), A.weights());
auto make_SqrtCauchyEigen_Te1(int exp, std::array< value_type, 2 > EVs, unsigned stepsCauchy)
Create a functor that computes using either Eigen or SqrtCauchy solve based on whichever is fastest ...
Definition: matrixfunction.h:131
This function uses an Eigen decomposition for small sizes of T and a SqrtCauchyInt solve for larger sizes to optimize execution times
- Parameters
-
exp | exponent if +1 compute \( \sqrt{T}\), if -1 compute \( 1/\sqrt{T}\) |
EVs | {minimum Eigenvalue of A, maximum Eigenvalue of A} |
stepsCauchy | iterations of cauchy integral |
- Returns
- an operator to use in
UniversalLanczos
solve method
- See also
SqrtCauchyInt
UniversalLanczos
◆ make_SqrtODE_Te1()
template<class value_type >
auto dg::mat::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.
- Parameters
-
exp | exponent if +1 compute \( \sqrt{T}\), if -1 compute \( 1/\sqrt{T}\) |
tableau | Tableau of time integrator |
rtol | relative tolerance of time integrator |
atol | absolute tolerance of time integrator |
number | links to number of steps in time integrator |
- Returns
- an operator to use in
UniversalLanczos
solve method
- See also
make_directODESolve
UniversalLanczos