Extension: Matrix functions
#include "dg/matrix/matrix.h"
Loading...
Searching...
No Matches
matrixfunction.h
Go to the documentation of this file.
1#pragma once
2#include <cmath>
3
4#include <boost/math/special_functions.hpp> // has to be included before lapack in certain versions
5#include "dg/algorithm.h"
6
7#include "functors.h"
8#include "sqrt_cauchy.h"
9#include "sqrt_ode.h"
10
11
12namespace dg {
13namespace mat {
16
28template<class UnaryOp>
29auto make_FuncEigen_Te1( UnaryOp f)
30{
31 return [f]( const auto& T)
32 {
33 using value_type = typename std::decay_t<decltype(T)>::value_type;
34 unsigned iter = T.O.size();
35 thrust::host_vector<value_type> evals = T.O, plus = T.P;
36 thrust::host_vector<value_type> work (2*iter-2);
38 //Compute Eigendecomposition
39 lapack::stev(LAPACK_COL_MAJOR, 'V', evals, plus, EHt.data(), work);
40 auto EH = EHt.transpose();
41 //for( unsigned u=0; u<iter; u++)
42 // std::cout << u << " "<<evals[u]<<std::endl;
43 //Compute f(T) e1 = E f(Lambda) E^t e1
44 dg::HVec e1H(iter,0.), yH(e1H);
45 e1H[0] = 1.;
46 yH.resize( iter);
47 dg::blas2::symv(EHt, e1H, yH);
48 dg::blas1::transform(evals, e1H, [f] (double x){
49 try{
50 return f(x);
51 }
52 catch(boost::exception& e) //catch boost overflow error
53 {
54 return 0.;
55 }
56 });
57 dg::blas1::pointwiseDot(e1H, yH, e1H);
58 dg::blas2::symv(EH, e1H, yH);
59 return yH;
60 };
61}
62
81template< class value_type>
82auto make_SqrtCauchy_Te1( int exp, std::array<value_type,2> EVs, unsigned stepsCauchy)
83{
84 return [=]( const auto& T)
85 {
86 unsigned size = T.O.size();
87 thrust::host_vector<value_type> e1H(size, 0.), yH(e1H);
88 e1H[0] = 1.;
89
90 dg::mat::SqrtCauchyInt<HVec> cauchysqrtH( e1H);
91 auto wTinv = [&w = cauchysqrtH.w(), &T = T]( const auto& y, auto& x)
92 {
93 // invert 1*T + w*wI
94 dg::mat::compute_Tinv_y( T, x, y, 1., w*w);
95 };
96 cauchysqrtH(T, wTinv, e1H, yH, EVs, stepsCauchy, exp);
97 return yH;
98 };
99}
100
120template< class value_type>
121auto make_SqrtCauchyEigen_Te1( int exp, std::array<value_type,2> EVs, unsigned stepsCauchy)
122{
123 std::function< value_type(value_type)> func = dg::SQRT<value_type>();
124 if( exp < 0)
125 func = [](value_type x){return 1./sqrt(x);};
126
127 auto eigen = make_FuncEigen_Te1( func);
128 auto cauchy = make_SqrtCauchy_Te1( exp, EVs, stepsCauchy);
129 return [=]( const auto& T)
130 {
131 unsigned size = T.O.size();
132 dg::HVec yH;
133 if ( size < 40)
134 yH = eigen( T);
135 else
136 yH = cauchy(T);
137 return yH;
138 };
139}
140
141
152template< class value_type>
153auto make_SqrtODE_Te1( int exp, std::string tableau, value_type rtol,
154 value_type atol, unsigned& number)
155{
156 return [=, &num = number](const auto& T)
157 {
158 unsigned size = T.O.size();
159 HVec e1H(size, 0), yH(e1H);
160 e1H[0] = 1.;
161 auto inv_sqrt = make_inv_sqrtodeTri( T, e1H);
162 auto sqrtHSolve = make_directODESolve( inv_sqrt,
163 tableau, rtol, atol, num);
164 dg::apply( sqrtHSolve, e1H, yH);
165 if( exp >= 0 )
166 {
167 dg::apply( T, yH, e1H);
168 return e1H;
169 }
170 return yH;
171 };
172}
173
180inline auto make_Linear_Te1( int exp)
181{
182 return [= ](const auto& T)
183 {
184 unsigned size = T.O.size();
185 HVec e1H(size, 0), yH(e1H);
186 e1H[0] = 1.;
187 if( exp < 0)
188 compute_Tinv_y( T, yH, e1H);
189 else
190 dg::blas2::symv( T, e1H, yH);
191 return yH;
192 };
193}
194
196} //namespace mat
197} //namespace dg
const std::vector< value_type > & data() const
SquareMatrix transpose() const
void pointwiseDot(value_type alpha, const ContainerType1 &x1, const ContainerType2 &x2, value_type1 beta, ContainerType &y)
void transform(const ContainerType1 &x, ContainerType &y, UnaryOp op)
void apply(get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y)
void symv(MatrixType &&M, const ContainerType1 &x, ContainerType2 &y)
auto make_directODESolve(ExplicitRHS &&ode, std::string tableau, value_type epsTimerel, value_type epsTimeabs, unsigned &number, value_type t0=0., value_type t1=1.)
Operator that integrates an ODE from 0 to 1 with an adaptive ERK class as timestepper
Definition sqrt_ode.h:38
void compute_Tinv_y(const dg::TriDiagonal< thrust::host_vector< value_type > > &T, thrust::host_vector< value_type > &x, const thrust::host_vector< value_type > &y, value_type a=1., value_type d=0.)
Computes the value of via Thomas algorithm.
Definition tridiaginv.h:124
auto make_FuncEigen_Te1(UnaryOp f)
Create a functor that uses Eigenvalue decomposition to compute for symmetric tridiagonal T.
Definition matrixfunction.h:29
auto make_SqrtCauchy_Te1(int exp, std::array< value_type, 2 > EVs, unsigned stepsCauchy)
Create a functor that computes using SqrtCauchyInt.
Definition matrixfunction.h:82
auto make_Linear_Te1(int exp)
Create a functor that computes directly.
Definition matrixfunction.h:180
auto make_SqrtODE_Te1(int exp, std::string tableau, value_type rtol, value_type atol, unsigned &number)
Create a functor that computes using ODE solve.
Definition matrixfunction.h:153
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:121
thrust::host_vector< double > HVec
InvSqrtODE< Container > make_inv_sqrtodeTri(const Matrix &TH, const Container &copyable)
Definition sqrt_ode.h:167
Classes for Krylov space approximations of a Matrix-Vector product.
Cauchy integral
Definition sqrt_cauchy.h:37
const double & w() const
The in .
Definition sqrt_cauchy.h:60