Fast computation of \( \vec x = A^{\pm 1/2}\vec b\) for self-adjoint positive definite \(A\).
More...
|
| | MatrixSqrt ()=default |
| | Construct empty.
|
| |
| template<class MatrixType > |
| | MatrixSqrt (MatrixType &A, int exp, const ContainerType &weights, value_type eps_rel, value_type nrmb_correction=1., unsigned max_iter=500, unsigned cauchy_steps=40) |
| | Construct from matrix.
|
| |
| template<class ... Params> |
| void | construct (Params &&...ps) |
| | Perfect forward parameters to one of the constructors.
|
| |
| unsigned | get_iter () const |
| | Get the number of Lanczos iterations in latest call to operator()
|
| |
| void | set_benchmark (bool benchmark, std::string message="SQRT") |
| | Set or unset performance timings during iterations.
|
| |
| template<class ContainerType0 , class ContainerType1 > |
| void | operator() (const ContainerType0 b, ContainerType1 &x) |
| | Apply matrix sqrt.
|
| |
template<class ContainerType>
struct dg::mat::MatrixSqrt< ContainerType >
Fast computation of \( \vec x = A^{\pm 1/2}\vec b\) for self-adjoint positive definite \(A\).
Convenience wrapper that uses dg::mat::UniversalLanczos combined with dg::mat::make_SqrtCauchyEigen_Te1 in its "universal" stopping criterion
- Note
- This is the fastest method to compute matrix square roots vector multiplications that we found to date
- Attention
- Just as in the Lanczos or PCG methods the matrix \( A\) needs to be positive-definite (i.e. it won't work for negative definite)
◆ container_type
template<class ContainerType >
◆ value_type
template<class ContainerType >
◆ MatrixSqrt() [1/2]
template<class ContainerType >
◆ MatrixSqrt() [2/2]
template<class ContainerType >
template<class MatrixType >
| dg::mat::MatrixSqrt< ContainerType >::MatrixSqrt |
( |
MatrixType & | A, |
|
|
int | exp, |
|
|
const ContainerType & | weights, |
|
|
value_type | eps_rel, |
|
|
value_type | nrmb_correction = 1., |
|
|
unsigned | max_iter = 500, |
|
|
unsigned | cauchy_steps = 40 ) |
|
inline |
Construct from matrix.
- Template Parameters
-
- Parameters
-
| A | self-adjoint matrix; is stored by reference |
| exp | exponent, if < 0 then 1/sqrt(A) is computed, else sqrt(A) |
| weights | the weights in which A is self-adjoint |
| eps_rel | relative accuracy of solution |
| nrmb_correction | absolute accuracy in units of eps_rel |
| max_iter | Maximum number of iterations in Lanczos tridiagonalization |
| cauchy_steps | number of cells in the Cauchy integral |
◆ construct()
template<class ContainerType >
template<class ... Params>
Perfect forward parameters to one of the constructors.
- Template Parameters
-
| Params | deduced by the compiler |
- Parameters
-
| ps | parameters forwarded to constructors |
◆ get_iter()
template<class ContainerType >
Get the number of Lanczos iterations in latest call to operator()
◆ operator()()
template<class ContainerType >
template<class ContainerType0 , class ContainerType1 >
| void dg::mat::MatrixSqrt< ContainerType >::operator() |
( |
const ContainerType0 | b, |
|
|
ContainerType1 & | x ) |
|
inline |
Apply matrix sqrt.
- Parameters
-
| b | input vector |
| x | output vector, contains \( x = A^{\pm 1/2} \vec b\) |
◆ set_benchmark()
template<class ContainerType >
| void dg::mat::MatrixSqrt< ContainerType >::set_benchmark |
( |
bool | benchmark, |
|
|
std::string | message = "SQRT" ) |
|
inline |
Set or unset performance timings during iterations.
- Parameters
-
| benchmark | If true, additional output will be written to std::cout during solution |
| message | An optional identifier that is printed together with the benchmark (intended use is to distinguish different messages in the output) |
The documentation for this struct was generated from the following file: