19 #ifndef LIBMESH_DENSE_MATRIX_IMPL_H
20 #define LIBMESH_DENSE_MATRIX_IMPL_H
27 #include "libmesh/dense_matrix.h"
28 #include "libmesh/dense_vector.h"
29 #include "libmesh/int_range.h"
30 #include "libmesh/libmesh.h"
32 #ifdef LIBMESH_HAVE_METAPHYSICL
35 template <
typename,
typename>
42 template <
typename T,
typename D>
44 template <
typename T,
typename D>
60 if (this->use_blas_lapack)
61 this->_multiply_blas(M2, LEFT_MULTIPLY);
74 this->resize (M2.
m(), M3.
n());
77 this->multiply(*
this, M2, M3);
97 this->resize (M2.
m(), M3.
n());
100 this->multiply(*
this, M2, M3);
108 if (this->use_blas_lapack)
109 this->_multiply_blas(
A, LEFT_MULTIPLY_TRANSPOSE);
123 const unsigned int n_rows =
A.m();
124 const unsigned int n_cols =
A.n();
127 this->resize(n_cols,n_cols);
130 for (
unsigned int i=0; i<n_cols; ++i)
131 for (
unsigned int j=0; j<=i; ++j)
132 for (
unsigned int k=0; k<n_rows; ++k)
133 (*
this)(i,j) +=
B(k,i)*
B(k,j);
136 for (
unsigned int i=0; i<n_cols; ++i)
137 for (
unsigned int j=i+1; j<n_cols; ++j)
138 (*
this)(i,j) = (*
this)(j,i);
145 this->resize (
A.n(),
B.n());
147 libmesh_assert_equal_to (
A.m(),
B.m());
148 libmesh_assert_equal_to (this->m(),
A.n());
149 libmesh_assert_equal_to (this->n(),
B.n());
151 const unsigned int m_s =
A.n();
152 const unsigned int p_s =
A.m();
153 const unsigned int n_s = this->n();
158 for (
unsigned int i=0; i<m_s; i++)
159 for (
unsigned int k=0; k<p_s; k++)
160 if (
A.transpose(i,k) != 0.)
161 for (
unsigned int j=0; j<n_s; j++)
162 (*
this)(i,j) +=
A.transpose(i,k)*
B(k,j);
171 template<
typename T2>
185 const unsigned int n_rows =
A.m();
186 const unsigned int n_cols =
A.n();
189 this->resize(n_cols,n_cols);
192 for (
unsigned int i=0; i<n_cols; ++i)
193 for (
unsigned int j=0; j<=i; ++j)
194 for (
unsigned int k=0; k<n_rows; ++k)
195 (*
this)(i,j) +=
B(k,i)*
B(k,j);
198 for (
unsigned int i=0; i<n_cols; ++i)
199 for (
unsigned int j=i+1; j<n_cols; ++j)
200 (*
this)(i,j) = (*
this)(j,i);
207 this->resize (
A.n(),
B.n());
209 libmesh_assert_equal_to (
A.m(),
B.m());
210 libmesh_assert_equal_to (this->m(),
A.n());
211 libmesh_assert_equal_to (this->n(),
B.n());
213 const unsigned int m_s =
A.n();
214 const unsigned int p_s =
A.m();
215 const unsigned int n_s = this->n();
220 for (
unsigned int i=0; i<m_s; i++)
221 for (
unsigned int k=0; k<p_s; k++)
222 if (
A.transpose(i,k) != 0.)
223 for (
unsigned int j=0; j<n_s; j++)
224 (*
this)(i,j) +=
A.transpose(i,k)*
B(k,j);
233 if (this->use_blas_lapack)
234 this->_multiply_blas(M3, RIGHT_MULTIPLY);
247 this->resize (M2.
m(), M3.
n());
249 this->multiply(*
this, M2, M3);
256 template<
typename T2>
269 this->resize (M2.
m(), M3.
n());
271 this->multiply(*
this, M2, M3);
280 if (this->use_blas_lapack)
281 this->_multiply_blas(
B, RIGHT_MULTIPLY_TRANSPOSE);
295 const unsigned int n_rows =
B.m();
296 const unsigned int n_cols =
B.n();
299 this->resize(n_rows,n_rows);
302 for (
unsigned int i=0; i<n_rows; ++i)
303 for (
unsigned int j=0; j<=i; ++j)
304 for (
unsigned int k=0; k<n_cols; ++k)
305 (*
this)(i,j) +=
A(i,k)*
A(j,k);
308 for (
unsigned int i=0; i<n_rows; ++i)
309 for (
unsigned int j=i+1; j<n_rows; ++j)
310 (*
this)(i,j) = (*
this)(j,i);
317 this->resize (
A.m(),
B.m());
319 libmesh_assert_equal_to (
A.n(),
B.n());
320 libmesh_assert_equal_to (this->m(),
A.m());
321 libmesh_assert_equal_to (this->n(),
B.m());
323 const unsigned int m_s =
A.m();
324 const unsigned int p_s =
A.n();
325 const unsigned int n_s = this->n();
330 for (
unsigned int j=0; j<n_s; j++)
331 for (
unsigned int k=0; k<p_s; k++)
332 if (
B.transpose(k,j) != 0.)
333 for (
unsigned int i=0; i<m_s; i++)
334 (*
this)(i,j) +=
A(i,k)*
B.transpose(k,j);
342 template<
typename T2>
356 const unsigned int n_rows =
B.m();
357 const unsigned int n_cols =
B.n();
360 this->resize(n_rows,n_rows);
363 for (
unsigned int i=0; i<n_rows; ++i)
364 for (
unsigned int j=0; j<=i; ++j)
365 for (
unsigned int k=0; k<n_cols; ++k)
366 (*
this)(i,j) +=
A(i,k)*
A(j,k);
369 for (
unsigned int i=0; i<n_rows; ++i)
370 for (
unsigned int j=i+1; j<n_rows; ++j)
371 (*
this)(i,j) = (*
this)(j,i);
378 this->resize (
A.m(),
B.m());
380 libmesh_assert_equal_to (
A.n(),
B.n());
381 libmesh_assert_equal_to (this->m(),
A.m());
382 libmesh_assert_equal_to (this->n(),
B.m());
384 const unsigned int m_s =
A.m();
385 const unsigned int p_s =
A.n();
386 const unsigned int n_s = this->n();
391 for (
unsigned int j=0; j<n_s; j++)
392 for (
unsigned int k=0; k<p_s; k++)
393 if (
B.transpose(k,j) != 0.)
394 for (
unsigned int i=0; i<m_s; i++)
395 (*
this)(i,j) +=
A(i,k)*
B.transpose(k,j);
407 libmesh_assert_equal_to (this->n(), arg.
size());
414 if(this->m() == 0 || this->n() == 0)
417 if (this->use_blas_lapack)
418 this->_matvec_blas(1., 0., dest, arg);
421 const unsigned int n_rows = this->m();
422 const unsigned int n_cols = this->n();
424 for (
unsigned int i=0; i<n_rows; i++)
425 for (
unsigned int j=0; j<n_cols; j++)
426 dest(i) += (*this)(i,j)*arg(j);
433 template<
typename T2>
438 libmesh_assert_equal_to (this->n(), arg.
size());
442 dest.resize(this->m());
445 if (this->m() == 0 || this->n() == 0)
448 const unsigned int n_rows = this->m();
449 const unsigned int n_cols = this->n();
451 for (
unsigned int i=0; i<n_rows; i++)
452 for (
unsigned int j=0; j<n_cols; j++)
453 dest(i) += (*this)(i,j)*arg(j);
463 libmesh_assert_equal_to (this->m(), arg.
size());
473 if (this->use_blas_lapack)
475 this->_matvec_blas(1., 0., dest, arg,
true);
479 const unsigned int n_rows = this->m();
480 const unsigned int n_cols = this->n();
488 for (
unsigned int i=0; i<n_cols; i++)
489 for (
unsigned int j=0; j<n_rows; j++)
490 dest(i) += (*this)(j,i)*arg(j);
497 template<
typename T2>
502 libmesh_assert_equal_to (this->m(), arg.
size());
506 dest.resize(this->n());
512 const unsigned int n_rows = this->m();
513 const unsigned int n_cols = this->n();
521 for (
unsigned int i=0; i<n_cols; i++)
522 for (
unsigned int j=0; j<n_rows; j++)
523 dest(i) += (*this)(j,i)*arg(j);
540 if (this->use_blas_lapack)
541 this->_matvec_blas(factor, 1., dest, arg);
545 this->vector_mult(temp, arg);
546 dest.
add(factor, temp);
553 template<
typename T2,
typename T3>
567 this->vector_mult(temp, arg);
568 dest.add(factor, temp);
580 dest.
resize(sub_m, sub_n);
581 for (
unsigned int i=0; i<sub_m; i++)
582 for (
unsigned int j=0; j<sub_n; j++)
583 dest(i,j) = (*this)(i,j);
592 const unsigned int m = a.
size();
593 const unsigned int n = b.
size();
596 for (
unsigned int i = 0; i < m; ++i)
597 for (
unsigned int j = 0; j < n; ++j)
606 get_principal_submatrix(sub_m, sub_m, dest);
614 dest.
resize(this->n(), this->m());
618 dest(i,j) = (*this)(j,i);
639 libmesh_assert_equal_to (this->m(), this->n());
641 switch(this->_decomposition_type)
645 if (this->use_blas_lapack)
646 this->_lu_decompose_lapack();
648 this->_lu_decompose ();
655 if (this->use_blas_lapack)
658 libmesh_fallthrough();
663 if (!(this->use_blas_lapack))
666 libmesh_fallthrough();
669 libmesh_error_msg(
"Error! This matrix already has a different decomposition...");
672 if (this->use_blas_lapack)
673 this->_lu_back_substitute_lapack (b, x);
675 this->_lu_back_substitute (b, x);
690 libmesh_assert_equal_to (this->m(), n_cols);
691 libmesh_assert_equal_to (this->m(), b.
size());
703 for (
unsigned int i=0; i<n_cols; ++i)
706 if (_pivots[i] != static_cast<pivot_index_t>(i))
711 for (
unsigned int j=0; j<i; ++j)
718 const unsigned int last_row = n_cols-1;
720 for (
int i=last_row; i>=0; --i)
722 for (
int j=i+1; j<static_cast<int>(n_cols); ++j)
739 libmesh_assert_equal_to (this->_decomposition_type, NONE);
748 _pivots.resize(n_rows);
750 for (
unsigned int i=0; i<n_rows; ++i)
757 for (
unsigned int j=i+1; j<n_rows; ++j)
760 if (the_max < candidate_max)
762 the_max = candidate_max;
772 if (_pivots[i] != static_cast<pivot_index_t>(i))
774 for (
unsigned int j=0; j<n_rows; ++j)
781 libmesh_error_msg(
"Matrix A is singular!");
785 const T diag_inv = 1. /
A(i,i);
786 for (
unsigned int j=i+1; j<n_rows; ++j)
800 for (
unsigned int row=i+1; row<n_rows; ++row)
801 for (
unsigned int col=i+1; col<n_rows; ++col)
802 A(row,col) -=
A(row,i) *
A(i,col);
807 this->_decomposition_type = LU;
826 _svd_lapack(sigma, U, VT);
836 _svd_solve_lapack(rhs, x, rcond);
846 _evd_lapack(lambda_real, lambda_imag);
857 _evd_lapack(lambda_real, lambda_imag, &VL,
nullptr);
868 _evd_lapack(lambda_real, lambda_imag,
nullptr, &VR);
880 _evd_lapack(lambda_real, lambda_imag, &VL, &VR);
888 switch(this->_decomposition_type)
895 if (this->use_blas_lapack)
896 this->_lu_decompose_lapack();
898 this->_lu_decompose ();
907 libmesh_error_msg(
"Error! Can't compute the determinant under the current decomposition.");
918 unsigned int n_interchanges = 0;
921 if (this->_decomposition_type==LU)
922 if (_pivots[i] != static_cast<pivot_index_t>(i))
926 if (this->_decomposition_type==LU_BLAS_LAPACK)
927 if (_pivots[i] != static_cast<pivot_index_t>(i+1))
930 determinant *= (*this)(i,i);
935 Real sign = n_interchanges % 2 == 0 ? 1. : -1.;
937 return sign*determinant;
945 template <
typename T>
946 template <
typename T2>
951 switch(this->_decomposition_type)
955 this->_cholesky_decompose ();
966 libmesh_error_msg(
"Error! This matrix already has a different decomposition...");
970 this->_cholesky_back_substitute (b, x);
983 libmesh_assert_equal_to (this->_decomposition_type, NONE);
991 libmesh_assert_equal_to (n_rows, n_cols);
996 for (
unsigned int i=0; i<n_rows; ++i)
998 for (
unsigned int j=i; j<n_cols; ++j)
1000 for (
unsigned int k=0; k<i; ++k)
1001 A(i,j) -=
A(i,k) *
A(j,k);
1005 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
1007 libmesh_error_msg(
"Error! Can only use Cholesky decomposition with symmetric positive definite matrices.");
1013 A(j,i) =
A(i,j) /
A(i,i);
1018 this->_decomposition_type = CHOLESKY;
1023 template <
typename T>
1024 template <
typename T2>
1034 libmesh_assert_equal_to (n_rows, n_cols);
1043 for (
unsigned int i=0; i<n_cols; ++i)
1047 for (
unsigned int k=0; k<i; ++k)
1048 temp -=
A(i,k)*x(k);
1050 x(i) = temp /
A(i,i);
1054 for (
unsigned int i=0; i<n_cols; ++i)
1056 const unsigned int ib = (n_cols-1)-i;
1058 for (
unsigned int k=(ib+1); k<n_cols; ++k)
1059 x(ib) -=
A(k,ib) * x(k);
1118 #endif // LIBMESH_DENSE_MATRIX_IMPL_H