Extension: Matrix functions
#include "dg/matrix/matrix.h"
dg::mat::MCGFuncEigen< Container > Class Template Reference

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...

Public Types

using value_type = dg::get_value_type< Container >
 
using HDiaMatrix = cusp::dia_matrix< int, value_type, cusp::host_memory >
 
using HCooMatrix = cusp::coo_matrix< int, value_type, cusp::host_memory >
 
using HArray2d = cusp::array2d< value_type, cusp::host_memory >
 
using HArray1d = cusp::array1d< value_type, cusp::host_memory >
 
using HVec = dg::HVec
 

Public Member Functions

 MCGFuncEigen ()
 Allocate nothing, Call construct method before usage. More...
 
 MCGFuncEigen (const Container &copyable, unsigned max_iterations)
 Construct MCGFuncEigen. More...
 
template<class ... Params>
void construct (Params &&...ps)
 Perfect forward parameters to one of the constructors. More...
 
template<class MatrixType , class ContainerType0 , class ContainerType1 , class ContainerType2 , class UnaryOp >
unsigned operator() (ContainerType0 &x, UnaryOp f, MatrixType &&A, const ContainerType1 &b, const ContainerType2 &weights, value_type eps, value_type nrmb_correction, value_type res_fac)
 Compute \(x \approx f(A)*b \) via M-CG method and Eigen decomposition. More...
 

Detailed Description

template<class Container>
class dg::mat::MCGFuncEigen< Container >

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.

Attention
This class has the exact same result as the corresponding UniversalLanczos class. Its use is for mere experimental purposes.

A is self-adjoint in the weights \( W\)

The approximation relies on Projection \(x = f(A) b \approx R f(\tilde T)^{-1} e_1\), where \(\tilde T\) and \(R\) are the tridiagonal and orthogonal matrix of the M-CG solve and \(e_1\) is the normalized unit vector. The vector \(f(\tilde T) e_1\) is computed via the Eigen decomposition and similarity transform of \( \tilde T\).

Note
Since MCG and Lanczos are equivalent the result of this class is the same within round-off errors as a Lanczos solve with the "residual" error norm

Member Typedef Documentation

◆ HArray1d

template<class Container >
using dg::mat::MCGFuncEigen< Container >::HArray1d = cusp::array1d< value_type, cusp::host_memory>

◆ HArray2d

template<class Container >
using dg::mat::MCGFuncEigen< Container >::HArray2d = cusp::array2d< value_type, cusp::host_memory>

◆ HCooMatrix

template<class Container >
using dg::mat::MCGFuncEigen< Container >::HCooMatrix = cusp::coo_matrix<int, value_type, cusp::host_memory>

◆ HDiaMatrix

template<class Container >
using dg::mat::MCGFuncEigen< Container >::HDiaMatrix = cusp::dia_matrix<int, value_type, cusp::host_memory>

◆ HVec

template<class Container >
using dg::mat::MCGFuncEigen< Container >::HVec = dg::HVec

◆ value_type

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

Constructor & Destructor Documentation

◆ MCGFuncEigen() [1/2]

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

Allocate nothing, Call construct method before usage.

◆ MCGFuncEigen() [2/2]

template<class Container >
dg::mat::MCGFuncEigen< Container >::MCGFuncEigen ( const Container &  copyable,
unsigned  max_iterations 
)
inline

Construct MCGFuncEigen.

Parameters
copyablea copyable container
max_iterationsMax iterations of Lanczos method

Member Function Documentation

◆ construct()

template<class Container >
template<class ... Params>
void dg::mat::MCGFuncEigen< 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 >
template<class MatrixType , class ContainerType0 , class ContainerType1 , class ContainerType2 , class UnaryOp >
unsigned dg::mat::MCGFuncEigen< Container >::operator() ( ContainerType0 &  x,
UnaryOp  f,
MatrixType &&  A,
const ContainerType1 &  b,
const ContainerType2 &  weights,
value_type  eps,
value_type  nrmb_correction,
value_type  res_fac 
)
inline

Compute \(x \approx f(A)*b \) via M-CG method and Eigen decomposition.

Parameters
binput vector
xoutput vector
fthe matrix function (e.g. dg::SQRT<double> or dg::EXP<double>)
Aself-adjoint and semi-positive definit matrix
weightsweights
epsrelative accuracy of M-CG method
nrmb_correctionthe absolute error C in units of eps to be respected
res_facresidual factor for stopping criterium of M-CG method
Returns
number of iterations of M-CG routine

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