22template<
class ContainerType>
43 template<
class MatrixType>
47 unsigned max_iter = 500,
unsigned cauchy_steps = 40
49 m_exp(exp), m_cauchy( cauchy_steps), m_eps(eps_rel),
50 m_abs(nrmb_correction)
52 m_A = [&](
const ContainerType&
x, ContainerType&
y){
55 m_lanczos.construct(
weights, max_iter);
61 template<
class ...Params>
65 *
this =
MatrixSqrt( std::forward<Params>( ps)...);
78 m_benchmark = benchmark;
88 template<
class ContainerType0,
class ContainerType1>
89 void operator()(
const ContainerType0 b, ContainerType1& x)
93 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
98 m_number = m_lanczos.solve(
x, func, m_A, b, m_weights, m_eps, m_abs,
102 DG_RANK0 std::cout <<
"# `"<<m_message<<
"` solve with {"<<m_number<<
","<<m_cauchy<<
"} iterations took "<<t.
diff()<<
"s\n";
106 ContainerType m_weights;
107 std::function< void(
const ContainerType&, ContainerType&)> m_A;
108 std::array<value_type, 2> m_EVs;
110 unsigned m_number, m_cauchy;
112 bool m_benchmark =
true;
113 std::string m_message =
"SQRT";
136template<
class ContainerType>
157 template<
class MatrixType>
161 unsigned max_iter = 500,
164 m_f_inner(f_inner), m_eps(eps_rel),
165 m_abs(nrmb_correction)
167 m_A = [&](
const ContainerType&
x, ContainerType&
y){
170 m_lanczos.construct(
weights, max_iter);
173 template<
class ...Params>
190 m_benchmark = benchmark;
201 template<
class UnaryOp,
class ContainerType0,
class ContainerType1>
202 void operator()( UnaryOp f_outer,
const ContainerType0 b, ContainerType1& x)
206 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
211 m_number = m_lanczos.solve(
x, func, m_A, b, m_weights, m_eps, m_abs,
215 DG_RANK0 std::cout <<
"# `"<<m_message<<
"` solve with {"<<m_number<<
"} iterations took "<<t.
diff()<<
"s\n";
219 ContainerType m_weights;
220 std::function< void(
const ContainerType&, ContainerType&)> m_A;
221 std::array<value_type, 2> m_EVs;
225 bool m_benchmark =
true;
226 std::string m_message =
"Function";
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.
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.