decompositions.hpp
Go to the documentation of this file.
1/*!
2 * \brief Defines types for matrix decompositions and solvers
3 *
4 * \author O. Krause
5 * \date 2016
6 *
7 *
8 * \par Copyright 1995-2015 Shark Development Team
9 *
10 * <BR><HR>
11 * This file is part of Shark.
12 * <http://image.diku.dk/shark/>
13 *
14 * Shark is free software: you can redistribute it and/or modify
15 * it under the terms of the GNU Lesser General Public License as published
16 * by the Free Software Foundation, either version 3 of the License, or
17 * (at your option) any later version.
18 *
19 * Shark is distributed in the hope that it will be useful,
20 * but WITHOUT ANY WARRANTY; without even the implied warranty of
21 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22 * GNU Lesser General Public License for more details.
23 *
24 * You should have received a copy of the GNU Lesser General Public License
25 * along with Shark. If not, see <http://www.gnu.org/licenses/>.
26 *
27 */
28#ifndef REMORA_DECOMPOSITIONS_HPP
29#define REMORA_DECOMPOSITIONS_HPP
30
31#include "kernels/trsm.hpp"
32#include "kernels/trsv.hpp"
33#include "kernels/potrf.hpp"
34#include "kernels/pstrf.hpp"
35#include "kernels/getrf.hpp"
36#include "kernels/syev.hpp"
37#include "assignment.hpp"
38#include "permutation.hpp"
39#include "matrix_expression.hpp"
40#include "proxy_expressions.hpp"
41#include "vector_expression.hpp"
42
43namespace remora{
44template<class D, class Device>
45struct solver_expression{
46 typedef Device device_type;
47
48 D const& operator()() const {
49 return *static_cast<D const*>(this);
50 }
51
52 D& operator()() {
53 return *static_cast<D*>(this);
54 }
55};
56
57
58/// \brief Lower triangular Cholesky decomposition.
59///
60/// Given an \f$ m \times m \f$ symmetric positive definite matrix
61/// \f$A\f$, represents the lower triangular matrix \f$L\f$ such that
62/// \f$A = LL^T \f$.
63///
64/// This decomposition is a corner block of many linear algebra routines
65/// especially for solving symmetric positive definite systems of equations.
66///
67/// Unlike many other decompositions, this decomposition is fast,
68/// numerically stable and can be updated when the original matrix is changed
69/// (rank-k updates of the form A<-alpha A + VCV^T)
70template<class MatrixStorage>
71class cholesky_decomposition:
72 public solver_expression<
73 cholesky_decomposition<MatrixStorage>,
74 typename MatrixStorage::device_type
75>{
76private:
77 typedef typename MatrixStorage::value_type value_type;
78public:
79 typedef typename MatrixStorage::device_type device_type;
80 cholesky_decomposition(){}
81 template<class E>
82 cholesky_decomposition(matrix_expression<E,device_type> const& e):m_cholesky(e){
83 kernels::potrf<lower>(m_cholesky);
84 }
85
86 template<class E>
87 void decompose(matrix_expression<E,device_type> const& e){
88 m_cholesky.resize(e().size1(), e().size2());
89 noalias(m_cholesky) = e;
90 kernels::potrf<lower>(m_cholesky);
91 }
92
93 MatrixStorage const& lower_factor()const{
94 return m_cholesky;
95 }
96
97 auto upper_factor()const -> decltype(trans(std::declval<MatrixStorage const&>())){
98 return trans(m_cholesky);
99 }
100
101 template<class MatB>
102 void solve(matrix_expression<MatB,device_type>& B, left)const{
103 kernels::trsm<lower,left >(lower_factor(),B);
104 kernels::trsm<upper,left >(upper_factor(),B);
105 }
106 template<class MatB>
107 void solve(matrix_expression<MatB,device_type>& B, right)const{
108 kernels::trsm<upper,right >(upper_factor(),B);
109 kernels::trsm<lower,right >(lower_factor(),B);
110 }
111 template<class VecB, bool Left>
112 void solve(vector_expression<VecB,device_type>&b, system_tag<Left>)const{
113 kernels::trsv<lower, left>(lower_factor(),b);
114 kernels::trsv<upper,left>(upper_factor(),b);
115 }
116
117 /// \brief Updates a covariance factor by a rank one update
118 ///
119 /// Let \f$ A=LL^T \f$ be a matrix with its lower cholesky factor. Assume we want to update
120 /// A using a simple rank-one update \f$ A = \alpha A+ \beta vv^T \f$. This invalidates L and
121 /// it needs to be recomputed which is O(n^3). instead we can update the factorisation
122 /// directly by performing a similar, albeit more complex algorithm on L, which can be done
123 /// in O(L^2).
124 ///
125 /// Alpha is not required to be positive, but if it is not negative, one has to be carefull
126 /// that the update would keep A positive definite. Otherwise the decomposition does not
127 /// exist anymore and an exception is thrown.
128 ///
129 /// \param v the update vector
130 /// \param alpha the scaling factor, must be positive.
131 /// \param beta the update factor. it Can be positive or negative
132 template<class VecV>
133 void update(value_type alpha, value_type beta, vector_expression<VecV,device_type> const& v){
134 if(beta == 0){
135 m_cholesky *= std::sqrt(alpha);
136 return;
137 }
138 //implementation blatantly stolen from Eigen
139 std::size_t n = v().size();
140 auto& L = m_cholesky;
141 typename vector_temporary<VecV>::type temp = v;
142 double beta_prime = 1;
143 double a = std::sqrt(alpha);
144 for(std::size_t j=0; j != n; ++j)
145 {
146 double Ljj = a * L(j,j);
147 double dj = Ljj * Ljj;
148 double wj = temp(j);
149 double swj2 = beta * wj * wj;
150 double gamma = dj * beta_prime + swj2;
151
152 double x = dj + swj2/beta_prime;
153 if (x <= 0.0)
154 throw std::invalid_argument("[cholesky_decomposition::update] update makes matrix indefinite, no update available");
155 double nLjj = std::sqrt(x);
156 L(j,j) = nLjj;
157 beta_prime += swj2/dj;
158
159 // Update the terms of L
160 if(j+1 <n){
161 subrange(column(L,j),j+1,n) *= a;
162 noalias(subrange(temp,j+1,n)) -= (wj/Ljj) * subrange(column(L,j),j+1,n);
163 if(gamma == 0)
164 continue;
165 subrange(column(L,j),j+1,n) *= nLjj/Ljj;
166 noalias(subrange(column(L,j),j+1,n))+= (nLjj * beta*wj/gamma)*subrange(temp,j+1,n);
167 }
168 }
169 }
170
171 template<typename Archive>
172 void serialize( Archive & ar, const std::size_t version ) {
173 ar & m_cholesky;
174 }
175private:
176 MatrixStorage m_cholesky;
177};
178
179
180/// \brief Symmetric eigenvalue decomposition A=QDQ^T
181///
182/// every symmetric matrix can be decomposed into its eigenvalue decomposition
183/// A=QDQ^T, where Q is an orthogonal matrix with Q^TQ=QQ^T=I
184/// and D is the diagonal matrix of eigenvalues of A.
185template<class MatrixStorage>
186class symm_eigenvalue_decomposition:
187 public solver_expression<
188 symm_eigenvalue_decomposition<MatrixStorage>,
189 typename MatrixStorage::device_type
190>{
191private:
192 typedef typename MatrixStorage::value_type value_type;
193 typedef typename vector_temporary<MatrixStorage>::type VectorStorage;
194public:
195 typedef typename MatrixStorage::device_type device_type;
196 symm_eigenvalue_decomposition(){}
197 template<class E>
198 symm_eigenvalue_decomposition(matrix_expression<E,device_type> const& e){
199 decompose(e);
200 }
201
202 template<class E>
203 void decompose(matrix_expression<E,device_type> const& e){
204 REMORA_SIZE_CHECK(e().size1() == e().size2());
205 m_eigenvectors.resize(e().size1(),e().size1());
206 m_eigenvalues.resize(e().size1());
207 noalias(m_eigenvectors) = e;
208
209 kernels::syev(m_eigenvectors,m_eigenvalues);
210 }
211
212 MatrixStorage const& Q()const{
213 return m_eigenvectors;
214 }
215 VectorStorage const& D()const{
216 return m_eigenvalues;
217 }
218
219
220 template<class MatB>
221 void solve(matrix_expression<MatB,device_type>& B, left)const{
222 B() = Q() % to_diagonal(elem_inv(D()))% trans(Q()) % B;
223 }
224 template<class MatB>
225 void solve(matrix_expression<MatB,device_type>& B, right)const{
226 auto transB = trans(B);
227 solve(transB,left());
228 }
229 template<class VecB>
230 void solve(vector_expression<VecB,device_type>&b, left)const{
231 b() = Q() % safe_div(trans(Q()) % b,D() ,0.0);
232 }
233
234 template<class VecB>
235 void solve(vector_expression<VecB,device_type>&b, right)const{
236 solve(b,left());
237 }
238
239 template<typename Archive>
240 void serialize( Archive & ar, const std::size_t version ) {
241 ar & m_eigenvectors;
242 ar & m_eigenvalues;
243 }
244private:
245 MatrixStorage m_eigenvectors;
246 VectorStorage m_eigenvalues;
247};
248
249
250
251
252template<class MatrixStorage>
253class pivoting_lu_decomposition:
254public solver_expression<
255 pivoting_lu_decomposition<MatrixStorage>,
256 typename MatrixStorage::device_type
257>{
258public:
259 typedef typename MatrixStorage::device_type device_type;
260 template<class E>
261 pivoting_lu_decomposition(matrix_expression<E,device_type> const& e)
262 :m_factor(e), m_permutation(e().size1()){
263 kernels::getrf(m_factor,m_permutation);
264 }
265
266 MatrixStorage const& factor()const{
267 return m_factor;
268 }
269
270 permutation_matrix const& permutation() const{
271 return m_permutation;
272 }
273
274 template<class MatB>
275 void solve(matrix_expression<MatB,device_type>& B, left)const{
276 swap_rows(m_permutation,B);
277 kernels::trsm<unit_lower,left>(m_factor,B);
278 kernels::trsm<upper,left>(m_factor,B);
279 }
280 template<class MatB>
281 void solve(matrix_expression<MatB,device_type>& B, right)const{
282 kernels::trsm<upper,right>(m_factor,B);
283 kernels::trsm<unit_lower,right>(m_factor,B);
284 swap_columns_inverted(m_permutation,B);
285 }
286 template<class VecB>
287 void solve(vector_expression<VecB,device_type>& b, left)const{
288 swap_rows(m_permutation,b);
289 kernels::trsv<unit_lower,left>(m_factor,b);
290 kernels::trsv<upper,left>(m_factor,b);
291 }
292 template<class VecB>
293 void solve(vector_expression<VecB,device_type>& b, right)const{
294 kernels::trsv<upper,right>(m_factor,b);
295 kernels::trsv<unit_lower,right>(m_factor,b);
296 swap_rows_inverted(m_permutation,b);
297 }
298private:
299 MatrixStorage m_factor;
300 permutation_matrix m_permutation;
301};
302
303
304// This is an implementation suggested by
305// "Fast Computation of Moore-Penrose Inverse Matrices"
306// applied to the special case of symmetric pos semi-def matrices
307// trading numerical accuracy vs speed. We go for speed.
308//
309// The fact that A is not full rank means it is not invertable,
310// so we solve it in a least squares sense.
311//
312// We use the formula for the pseudo-inverse:
313// (P^T A P)^-1 = L(L^TL)^-1(L^TL)^-1 L^T
314// where L is a matrix obtained by some rank revealing factorization
315// P^T A P = L L^T
316// we chose a pivoting cholesky to make use of the fact that A is symmetric
317// and all eigenvalues are >=0. If A has full rank, this reduces to
318// the cholesky factor where the pivoting leads to slightly smaller numerical errors
319// At a higher computational cost compared to the normal cholesky decomposition.
320template<class MatrixStorage>
321class symm_pos_semi_definite_solver:
322 public solver_expression<symm_pos_semi_definite_solver<MatrixStorage>, typename MatrixStorage::device_type>{
323public:
324 typedef typename MatrixStorage::device_type device_type;
325 template<class E>
326 symm_pos_semi_definite_solver(matrix_expression<E,device_type> const& e)
327 :m_factor(e), m_permutation(e().size1()){
328 m_rank = kernels::pstrf<lower>(m_factor,m_permutation);
329 if(m_rank == e().size1()) return; //full rank, so m_factor is lower triangular and we are done
330
331 auto L = columns(m_factor,0,m_rank);
332 m_cholesky.decompose(prod(trans(L),L));
333 }
334
335 std::size_t rank()const{
336 return m_rank;
337 }
338
339 //compute C so that A^dagger = CC^T
340 //where A^dagger is the moore-penrose inverse
341 // m must be of size rank x n
342 template<class Mat>
343 void compute_inverse_factor(matrix_expression<Mat,device_type>& C)const{
344 REMORA_SIZE_CHECK(C().size1() == m_rank);
345 REMORA_SIZE_CHECK(C().size2() == m_factor.size1());
346 if(m_rank == m_factor.size1()){//matrix has full rank
347 //initialize as identity matrix and solve
348 noalias(C) = identity_matrix<double>( m_factor.size1());
349 swap_columns_inverted(m_permutation,C);
350 kernels::trsm<lower,left>(m_factor,C);
351 }else{
352 auto L = columns(m_factor,0,m_rank);
353 noalias(C) = trans(L);
354 m_cholesky.solve(C,left());
355 swap_columns_inverted(m_permutation,C);
356 }
357
358 }
359 template<class MatB>
360 void solve(matrix_expression<MatB,device_type>& B, left)const{
361 swap_rows(m_permutation,B);
362 if(m_rank == 0){//matrix is zero
363 B().clear();
364 }else if(m_rank == m_factor.size1()){//matrix has full rank
365 kernels::trsm<lower,left>(m_factor,B);
366 kernels::trsm<upper,left>(trans(m_factor),B);
367 }else{//matrix is missing rank
368 auto L = columns(m_factor,0,m_rank);
369 auto Z = eval_block(prod(trans(L),B));
370 m_cholesky.solve(Z,left());
371 m_cholesky.solve(Z,left());
372 noalias(B) = prod(L,Z);
373 }
374 swap_rows_inverted(m_permutation,B);
375 }
376 template<class MatB>
377 void solve(matrix_expression<MatB,device_type>& B, right)const{
378 //compute using symmetry of the system of equations
379 auto transB = trans(B);
380 solve(transB, left());
381 }
382 template<class VecB, bool Left>
383 void solve(vector_expression<VecB,device_type>& b, system_tag<Left>)const{
384 swap_rows(m_permutation,b);
385 if(m_rank == 0){//matrix is zero
386 b().clear();
387 }else if(m_rank == m_factor.size1()){//matrix has full rank
388 kernels::trsv<lower,left >(m_factor,b);
389 kernels::trsv<upper,left >(trans(m_factor),b);
390 }else{//matrix is missing rank
391 auto L = columns(m_factor,0,m_rank);
392 auto z = eval_block(prod(trans(L),b));
393 m_cholesky.solve(z,left());
394 m_cholesky.solve(z,left());
395 noalias(b) = prod(L,z);
396 }
397 swap_rows_inverted(m_permutation,b);
398 }
399private:
400 std::size_t m_rank;
401 MatrixStorage m_factor;
402 cholesky_decomposition<MatrixStorage> m_cholesky;
403 permutation_matrix m_permutation;
404};
405
406template<class MatA>
407class cg_solver:public solver_expression<cg_solver<MatA>, typename MatA::device_type>{
408public:
409 typedef typename MatA::const_closure_type matrix_closure_type;
410 typedef typename MatA::device_type device_type;
411 cg_solver(matrix_closure_type const& e, double epsilon = 1.e-10, unsigned int max_iterations = 0)
412 :m_expression(e), m_epsilon(epsilon), m_max_iterations(max_iterations){}
413
414 template<class MatB>
415 void solve(
416 matrix_expression<MatB,device_type>& B, left,
417 double epsilon, unsigned max_iterations
418 )const{
419 typename matrix_temporary<MatB>::type X = B;
420 cg(m_expression,X, B, epsilon, max_iterations);
421 noalias(B) = X;
422 }
423 template<class MatB>
424 void solve(
425 matrix_expression<MatB,device_type>& B, right,
426 double epsilon, unsigned max_iterations
427 )const{
428 auto transB = trans(B);
429 typename transposed_matrix_temporary<MatB>::type X = transB;
430 cg(m_expression,X,transB, epsilon, max_iterations);
431 noalias(transB) = X;
432 }
433 template<class VecB, bool Left>
434 void solve(
435 vector_expression<VecB,device_type>&b, system_tag<Left>,
436 double epsilon, unsigned max_iterations
437 )const{
438 typename vector_temporary<VecB>::type x = b;
439 cg(m_expression,x,b,epsilon,max_iterations);
440 noalias(b) = x;
441 }
442
443 template<class VecB, bool Left>
444 void solve(vector_expression<VecB,device_type>&b, system_tag<Left> tag)const{
445 solve(b, tag, m_epsilon, m_max_iterations);
446 }
447 template<class MatB, bool Left>
448 void solve(matrix_expression<MatB,device_type>&B, system_tag<Left> tag)const{
449 solve(B, tag, m_epsilon, m_max_iterations);
450 }
451private:
452 template<class MatT, class VecB, class VecX>
453 void cg(
454 matrix_expression<MatT, device_type> const& A,
455 vector_expression<VecX, device_type>& x,
456 vector_expression<VecB, device_type> const& b,
457 double epsilon,
458 unsigned int max_iterations
459 )const{
460 REMORA_SIZE_CHECK(A().size1() == A().size2());
461 REMORA_SIZE_CHECK(A().size1() == b().size());
462 REMORA_SIZE_CHECK(A().size1() == x().size());
463 typedef typename vector_temporary<VecX>::type vector_type;
464
465 std::size_t dim = b().size();
466
467 //initialize point.
468 vector_type residual = b - prod(A,x);
469 //check if provided solution is better than starting at 0
470 if(norm_inf(residual) > norm_inf(b)){
471 x().clear();
472 residual = b;
473 }
474 vector_type next_residual(dim); //the next residual
475 vector_type p = residual; //the search direction- initially it is the gradient direction
476 vector_type Ap(dim); //stores prod(A,p)
477
478 for(std::size_t iter = 0;; ++iter){
479 if(max_iterations != 0 && iter >= max_iterations) break;
480 noalias(Ap) = prod(A,p);
481 double rsqr = norm_sqr(residual);
482 double alpha = rsqr/inner_prod(p,Ap);
483 noalias(x) += alpha * p;
484 noalias(next_residual) = residual - alpha * Ap;
485 if(norm_inf(next_residual) < epsilon)
486 break;
487
488 double beta = inner_prod(next_residual,next_residual)/rsqr;
489 p *= beta;
490 noalias(p) += next_residual;
491 swap(residual,next_residual);
492 }
493 }
494
495 template<class MatT, class MatB, class MatX>
496 void cg(
497 matrix_expression<MatT, device_type> const& A,
498 matrix_expression<MatX, device_type>& X,
499 matrix_expression<MatB, device_type> const& B,
500 double epsilon,
501 unsigned int max_iterations
502 )const{
503 REMORA_SIZE_CHECK(A().size1() == A().size2());
504 REMORA_SIZE_CHECK(A().size1() == B().size1());
505 REMORA_SIZE_CHECK(A().size1() == X().size1());
506 REMORA_SIZE_CHECK(B().size2() == X().size2());
507 typedef typename vector_temporary<MatX>::type vector_type;
508 typedef typename matrix_temporary<MatX>::type matrix_type;
509
510 std::size_t dim = B().size1();
511 std::size_t num_rhs = B().size2();
512
513 //initialize gradient given the starting point
514 matrix_type residual = B - prod(A,X);
515 //check for each rhs whether the starting point is better than starting from scratch
516 for(std::size_t i = 0; i != num_rhs; ++i){
517 if(norm_inf(column(residual,i)) <= norm_inf(column(residual,i))){
518 column(X,i).clear();
519 noalias(column(residual,i)) = column(B,i);
520 }
521 }
522
523 vector_type next_residual(dim); //the next residual of a column
524 matrix_type P = residual; //the search direction- initially it is the gradient direction
525 matrix_type AP(dim, num_rhs); //stores prod(A,p)
526
527 for(std::size_t iter = 0;; ++iter){
528 if(max_iterations != 0 && iter >= max_iterations) break;
529 //compute the product for all rhs at the same time
530 noalias(AP) = prod(A,P);
531 //for each rhs apply a step of cg
532 for(std::size_t i = 0; i != num_rhs; ++i){
533 auto r = column(residual,i);
534 //skip this if we are done already
535 //otherwise we might run into numerical troubles later on
536 if(norm_inf(r) < epsilon) continue;
537
538 auto x = column(X,i);
539 auto p = column(P,i);
540 auto Ap = column(AP,i);
541 double rsqr = norm_sqr(r);
542 double alpha = rsqr/inner_prod(p,Ap);
543 noalias(x) += alpha * p;
544 noalias(next_residual) = r - alpha * Ap;
545 double beta = inner_prod(next_residual,next_residual)/rsqr;
546 p *= beta;
547 noalias(p) += next_residual;
548 noalias(r) = next_residual;
549 }
550 //if all solutions are within tolerance, we are done
551 if(max(abs(residual)) < epsilon)
552 break;
553 }
554 }
555
556 matrix_closure_type m_expression;
557 double m_epsilon;
558 unsigned int m_max_iterations;
559};
560
561
562/////////////////////////////////////////////////////////////////
563////Traits connecting decompositions/solvers with tags
564/////////////////////////////////////////////////////////////////
565struct symm_pos_def{ typedef symm_pos_def transposed_orientation;};
566struct symm_semi_pos_def{ typedef symm_semi_pos_def transposed_orientation;};
567struct indefinite_full_rank{ typedef indefinite_full_rank transposed_orientation;};
568struct conjugate_gradient{
569 typedef conjugate_gradient transposed_orientation;
570 double epsilon;
571 unsigned max_iterations;
572 conjugate_gradient(double epsilon = 1.e-10, unsigned max_iterations = 0)
573 :epsilon(epsilon), max_iterations(max_iterations){}
574};
575
576
577namespace detail{
578template<class MatA, class SolverTag>
579struct solver_traits;
580
581template<class MatA, bool Upper, bool Unit>
582struct solver_traits<MatA,triangular_tag<Upper,Unit> >{
583 class type: public solver_expression<type,typename MatA::device_type>{
584 public:
585 typedef typename MatA::const_closure_type matrix_closure_type;
586 typedef typename MatA::device_type device_type;
587 type(matrix_closure_type const& e, triangular_tag<Upper,Unit>):m_matrix(e){}
588
589 template<class MatB, bool Left>
590 void solve(matrix_expression<MatB, device_type>& B, system_tag<Left> ){
591 kernels::trsm<triangular_tag<Upper,Unit>,system_tag<Left> >(m_matrix,B);
592 }
593 template<class VecB, bool Left>
594 void solve(vector_expression<VecB, device_type>& b, system_tag<Left> ){
595 kernels::trsv<triangular_tag<Upper,Unit>,system_tag<Left> >(m_matrix,b);
596 }
597 private:
598 matrix_closure_type m_matrix;
599 };
600};
601
602template<class MatA>
603struct solver_traits<MatA,symm_pos_def>{
604 struct type : public cholesky_decomposition<typename matrix_temporary<MatA>::type>{
605 template<class M>
606 type(M const& m, symm_pos_def)
607 :cholesky_decomposition<typename matrix_temporary<MatA>::type>(m){}
608 };
609};
610
611template<class MatA>
612struct solver_traits<MatA,indefinite_full_rank>{
613 struct type : public pivoting_lu_decomposition<typename matrix_temporary<MatA>::type>{
614 template<class M>
615 type(M const& m, indefinite_full_rank)
616 :pivoting_lu_decomposition<typename matrix_temporary<MatA>::type>(m){}
617 };
618};
619
620template<class MatA>
621struct solver_traits<MatA,symm_semi_pos_def>{
622 struct type : public symm_pos_semi_definite_solver<typename matrix_temporary<MatA>::type>{
623 template<class M>
624 type(M const& m, symm_semi_pos_def)
625 :symm_pos_semi_definite_solver<typename matrix_temporary<MatA>::type>(m){}
626 };
627};
628
629template<class MatA>
630struct solver_traits<MatA,conjugate_gradient>{
631 struct type : public cg_solver<MatA>{
632 template<class M>
633 type(M const& m, conjugate_gradient t):cg_solver<MatA>(m,t.epsilon,t.max_iterations){}
634 };
635};
636
637}
638
639template<class MatA, class SolverType>
640struct solver:public detail::solver_traits<MatA,SolverType>::type{
641 template<class M>
642 solver(M const& m, SolverType t = SolverType()): detail::solver_traits<MatA,SolverType>::type(m,t){}
643};
644
645
646}
647#endif