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

Right hand side of the square root ODE. More...

Public Types

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

Public Member Functions

 InvSqrtODE ()
 
template<class MatrixType >
 InvSqrtODE (MatrixType &A, const Container &copyable)
 Construct SqrtOde operator. More...
 
template<class ... Params>
void construct (Params &&...ps)
 Perfect forward parameters to one of the constructors. More...
 
const value_typetime () const
 
auto make_operator () const
 
template<class MatrixType >
void set_inverse_operator (const MatrixType &OpInv)
 
void operator() (value_type t, const Container &y, Container &yp)
 Compute rhs term. More...
 

Detailed Description

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

Right hand side of the square root ODE.

\[ \dot{y}= \left[(1-t) A +t I\right]^{-1} (I - A)/2 y \]

where \( A\) is the matrix

This class is based on the approach of the paper Numerical approximation of the product of the square root of a matrix with a vector by E. J. Allen et al

Note
Solution of ODE: \( y(1) = \frac{1}{\sqrt{A}} y(0)\) If \( y(0) = A b \) then we get \( y(1) = \sqrt{A} b\)

Member Typedef Documentation

◆ container_type

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

◆ value_type

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

Constructor & Destructor Documentation

◆ InvSqrtODE() [1/2]

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

◆ InvSqrtODE() [2/2]

template<class Container >
template<class MatrixType >
dg::mat::InvSqrtODE< Container >::InvSqrtODE ( MatrixType &  A,
const Container &  copyable 
)
inline

Construct SqrtOde operator.

Parameters
Amatrix (stored by reference so needs to live)
copyablecopyable container

Member Function Documentation

◆ construct()

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

Perfect forward parameters to one of the constructors.

Template Parameters
Paramsdeduced by the compiler
Parameters
psparameters forwarded to constructors

◆ make_operator()

template<class Container >
auto dg::mat::InvSqrtODE< Container >::make_operator ( ) const
inline

◆ operator()()

template<class Container >
void dg::mat::InvSqrtODE< Container >::operator() ( value_type  t,
const Container &  y,
Container &  yp 
)
inline

Compute rhs term.

i.e.

\[ yp= (tI + (1-t) A)^{-1} (I - A)/2 y \]

Parameters
tis time
yis \( y\)
ypis \( \dot{y}\)
Note
Solution of ODE: \( y(1) = 1/\sqrt{A} y(0)\)

◆ set_inverse_operator()

template<class Container >
template<class MatrixType >
void dg::mat::InvSqrtODE< Container >::set_inverse_operator ( const MatrixType &  OpInv)
inline

◆ time()

template<class Container >
const value_type & dg::mat::InvSqrtODE< Container >::time ( ) const
inline

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