Extension: Matrix functions
#include "dg/matrix/matrix.h"
Loading...
Searching...
No Matches
dg::mat::DirectSqrtCauchy< Container > Struct Template Reference

Shortcut for \(b \approx \sqrt{A} x \) solve directly via sqrt Cauchy combined with PCG inversions. More...

Public Types

using container_type = Container
 
using value_type = dg::get_value_type<Container>
 

Public Member Functions

 DirectSqrtCauchy ()
 empty object ( no memory allocation)
 
template<class MatrixType >
 DirectSqrtCauchy (MatrixType &A, const Container &weights, value_type epsCG, unsigned iterCauchy, std::array< value_type, 2 > EVs, int exp)
 Construct DirectSqrtCauchy.
 
template<class ... Params>
void construct (Params &&...ps)
 Perfect forward parameters to one of the constructors.
 
unsigned operator() (const Container &b, Container &x)
 Compute \(x \approx \sqrt{A} b \) via sqrt Cauchy integral solve.
 

Detailed Description

template<class Container>
struct dg::mat::DirectSqrtCauchy< Container >

Shortcut for \(b \approx \sqrt{A} x \) solve directly via sqrt Cauchy combined with PCG inversions.

Member Typedef Documentation

◆ container_type

template<class Container >
using dg::mat::DirectSqrtCauchy< Container >::container_type = Container

◆ value_type

template<class Container >
using dg::mat::DirectSqrtCauchy< Container >::value_type = dg::get_value_type<Container>

Constructor & Destructor Documentation

◆ DirectSqrtCauchy() [1/2]

template<class Container >
dg::mat::DirectSqrtCauchy< Container >::DirectSqrtCauchy ( )
inline

empty object ( no memory allocation)

◆ DirectSqrtCauchy() [2/2]

template<class Container >
template<class MatrixType >
dg::mat::DirectSqrtCauchy< Container >::DirectSqrtCauchy ( MatrixType & A,
const Container & weights,
value_type epsCG,
unsigned iterCauchy,
std::array< value_type, 2 > EVs,
int exp )
inline

Construct DirectSqrtCauchy.

Note
The Eigenvalues can be estimated from a few lanczos iterations (which is at least more reliable than doing it semi-analytically)
dg::mat::UniversalLanczos lanczos( A.weights(), 20);
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:66
std::array< value_type, 2 > compute_extreme_EV(const dg::TriDiagonal< thrust::host_vector< value_type > > &T)
Compute extreme Eigenvalues of a symmetric tridiangular matrix.
Definition tridiaginv.h:665
Parameters
AThe matrix (stored by reference)
weights
epsCG
iterCauchymaximum number of Cauchy iterations
EVs{minimum Eigenvalue of A, maximum Eigenvalue of A}
expif < 0 then the inverse sqrt is computed, else the sqrt

Member Function Documentation

◆ construct()

template<class Container >
template<class ... Params>
void dg::mat::DirectSqrtCauchy< Container >::construct ( Params &&... ps)
inline

Perfect forward parameters to one of the constructors.

Template Parameters
Paramsdeduced by the compiler
Parameters
psparameters forwarded to constructors

◆ operator()()

template<class Container >
unsigned dg::mat::DirectSqrtCauchy< Container >::operator() ( const Container & b,
Container & x )
inline

Compute \(x \approx \sqrt{A} b \) via sqrt Cauchy integral solve.

Parameters
binput vector
xoutput vector. Is approximating \(x \approx \sqrt{A} b \)
Returns
number of integration steps of sqrt cauchy solve

The documentation for this struct was generated from the following file: