Extension: Matrix functions
#include "dg/matrix/matrix.h"
matrixsqrt.h
Go to the documentation of this file.
1#pragma once
2
3#include "dg/algorithm.h"
4#include "lanczos.h"
5
6namespace dg{
7namespace mat{
8
9
22template<class ContainerType>
24{
25 using container_type = ContainerType;
27
29 MatrixSqrt() = default;
30
43 template<class MatrixType>
44 MatrixSqrt( MatrixType& A, int exp,
45 const ContainerType& weights, value_type eps_rel,
46 value_type nrmb_correction = 1.,
47 unsigned max_iter = 500, unsigned cauchy_steps = 40
48 ) : m_weights(weights),
49 m_exp(exp), m_cauchy( cauchy_steps), m_eps(eps_rel),
50 m_abs(nrmb_correction)
51 {
52 m_A = [&]( const ContainerType& x, ContainerType& y){
53 return dg::apply( A, x, y);
54 };
55 m_lanczos.construct( weights, max_iter);
57 auto T = eigen.tridiag( A, weights, weights);
59 }
61 template<class ...Params>
62 void construct( Params&& ...ps)
63 {
64 //construct and swap
65 *this = MatrixSqrt( std::forward<Params>( ps)...);
66 }
68 unsigned get_iter() const{return m_number;}
69
77 void set_benchmark( bool benchmark, std::string message = "SQRT"){
78 m_benchmark = benchmark;
79 m_message = message;
80 }
81
88 template<class ContainerType0, class ContainerType1>
89 void operator()( const ContainerType0 b, ContainerType1& x)
90 {
91#ifdef MPI_VERSION
92 int rank;
93 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
94#endif //MPI
95 dg::Timer t;
96 t.tic();
97 auto func = make_SqrtCauchyEigen_Te1( m_exp, m_EVs, m_cauchy);
98 m_number = m_lanczos.solve( x, func, m_A, b, m_weights, m_eps, m_abs,
99 "universal", 1., 2);
100 t.toc();
101 if( m_benchmark)
102 DG_RANK0 std::cout << "# `"<<m_message<<"` solve with {"<<m_number<<","<<m_cauchy<<"} iterations took "<<t.diff()<<"s\n";
103 }
104 private:
106 ContainerType m_weights;
107 std::function< void( const ContainerType&, ContainerType&)> m_A;
108 std::array<value_type, 2> m_EVs;
109 int m_exp;
110 unsigned m_number, m_cauchy;
111 value_type m_eps, m_abs;
112 bool m_benchmark = true;
113 std::string m_message = "SQRT";
114
115};
116
117// The following is only indirectly tested in diffusion project but should appear formally in _t file here as well
136template<class ContainerType>
138{
139 using container_type = ContainerType;
141
143 MatrixFunction() = default;
144
157 template<class MatrixType>
158 MatrixFunction( MatrixType& A,
159 const ContainerType& weights, value_type eps_rel,
160 value_type nrmb_correction = 1.,
161 unsigned max_iter = 500,
162 std::function<value_type(value_type)> f_inner = [](value_type x){return x;}
163 ) : m_weights(weights),
164 m_f_inner(f_inner), m_eps(eps_rel),
165 m_abs(nrmb_correction)
166 {
167 m_A = [&]( const ContainerType& x, ContainerType& y){
168 return dg::apply( A, x, y);
169 };
170 m_lanczos.construct( weights, max_iter);
171 }
173 template<class ...Params>
174 void construct( Params&& ...ps)
175 {
176 //construct and swap
177 *this = MatrixFunction( std::forward<Params>( ps)...);
178 }
180 unsigned get_iter() const{return m_number;}
181
189 void set_benchmark( bool benchmark, std::string message = "Function"){
190 m_benchmark = benchmark;
191 m_message = message;
192 }
193
201 template<class UnaryOp, class ContainerType0, class ContainerType1>
202 void operator()( UnaryOp f_outer, const ContainerType0 b, ContainerType1& x)
203 {
204#ifdef MPI_VERSION
205 int rank;
206 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
207#endif //MPI
208 dg::Timer t;
209 t.tic();
210 auto func = make_FuncEigen_Te1( [&](value_type x) {return f_outer(m_f_inner(x));});
211 m_number = m_lanczos.solve( x, func, m_A, b, m_weights, m_eps, m_abs,
212 "universal", 1., 2);
213 t.toc();
214 if( m_benchmark)
215 DG_RANK0 std::cout << "# `"<<m_message<<"` solve with {"<<m_number<<"} iterations took "<<t.diff()<<"s\n";
216 }
217 private:
219 ContainerType m_weights;
220 std::function< void( const ContainerType&, ContainerType&)> m_A;
221 std::array<value_type, 2> m_EVs;
222 std::function<value_type(value_type)> m_f_inner;
223 unsigned m_number;
224 value_type m_eps, m_abs;
225 bool m_benchmark = true;
226 std::string m_message = "Function";
227
228};
229}//namespace mat
230}//namespace dg
Tridiagonalize and approximate via Lanczos algorithm. A is self-adjoint in the weights .
Definition: lanczos.h:68
const HDiaMatrix & tridiag(MatrixType &&A, const ContainerType0 &b, const ContainerType1 &weights, value_type eps=1e-12, value_type nrmb_correction=1., std::string error_norm="universal", value_type res_fac=1., unsigned q=1)
Tridiagonalization of A using Lanczos method.
Definition: lanczos.h:184
void apply(get_value_type< ContainerType1 > alpha, MatrixType &&M, const ContainerType1 &x, get_value_type< ContainerType1 > beta, ContainerType2 &y)
MPI_Vector< thrust::host_vector< real_type > > weights(const aRealMPITopology2d< real_type > &g)
std::array< value_type, 2 > compute_extreme_EV(const cusp::dia_matrix< int, value_type, cusp::host_memory > &T)
Compute extreme Eigenvalues of a symmetric tridiangular matrix.
Definition: tridiaginv.h:631
typename TensorTraits< std::decay_t< Vector > >::value_type get_value_type
auto make_FuncEigen_Te1(UnaryOp f)
Create a functor that uses Eigenvalue decomposition to compute for symmetric tridiagonal T.
Definition: matrixfunction.h:35
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:131
Classes for Krylov space approximations of a Matrix-Vector product.
double diff() const
void toc()
void tic()
Computation of for self-adjoint positive definite .
Definition: matrixsqrt.h:138
void operator()(UnaryOp f_outer, const ContainerType0 b, ContainerType1 &x)
Apply matrix function.
Definition: matrixsqrt.h:202
MatrixFunction(MatrixType &A, const ContainerType &weights, value_type eps_rel, value_type nrmb_correction=1., unsigned max_iter=500, std::function< value_type(value_type)> f_inner=[](value_type x){return x;})
Construct from matrix.
Definition: matrixsqrt.h:158
unsigned get_iter() const
Get the number of Lanczos iterations in latest call to operator()
Definition: matrixsqrt.h:180
dg::get_value_type< ContainerType > value_type
Definition: matrixsqrt.h:140
ContainerType container_type
Definition: matrixsqrt.h:139
MatrixFunction()=default
Construct empty.
void set_benchmark(bool benchmark, std::string message="Function")
Set or unset performance timings during iterations.
Definition: matrixsqrt.h:189
void construct(Params &&...ps)
Perfect forward parameters to one of the constructors.
Definition: matrixsqrt.h:174
Fast computation of for self-adjoint positive definite .
Definition: matrixsqrt.h:24
void operator()(const ContainerType0 b, ContainerType1 &x)
Apply matrix sqrt.
Definition: matrixsqrt.h:89
void construct(Params &&...ps)
Perfect forward parameters to one of the constructors.
Definition: matrixsqrt.h:62
dg::get_value_type< ContainerType > value_type
Definition: matrixsqrt.h:26
unsigned get_iter() const
Get the number of Lanczos iterations in latest call to operator()
Definition: matrixsqrt.h:68
ContainerType container_type
Definition: matrixsqrt.h:25
void set_benchmark(bool benchmark, std::string message="SQRT")
Set or unset performance timings during iterations.
Definition: matrixsqrt.h:77
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.
Definition: matrixsqrt.h:44
MatrixSqrt()=default
Construct empty.
#define DG_RANK0