dense_matrix.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_DENSE_MATRIX_H
21 #define LIBMESH_DENSE_MATRIX_H
22 
23 // Local Includes
24 #include "libmesh/libmesh_common.h"
26 
27 // For the definition of PetscBLASInt.
28 #if (LIBMESH_HAVE_PETSC)
29 # include "libmesh/petsc_macro.h"
30 # include <petscsys.h>
31 #endif
32 
33 // C++ includes
34 #include <vector>
35 #include <algorithm>
36 
37 namespace libMesh
38 {
39 
40 // Forward Declarations
41 template <typename T> class DenseVector;
42 
53 template<typename T>
54 class DenseMatrix : public DenseMatrixBase<T>
55 {
56 public:
57 
61  DenseMatrix(const unsigned int new_m=0,
62  const unsigned int new_n=0);
63 
68  DenseMatrix (DenseMatrix &&) = default;
69  DenseMatrix (const DenseMatrix &) = default;
70  DenseMatrix & operator= (const DenseMatrix &) = default;
71  DenseMatrix & operator= (DenseMatrix &&) = default;
72  virtual ~DenseMatrix() = default;
73 
74  virtual void zero() override;
75 
79  T operator() (const unsigned int i,
80  const unsigned int j) const;
81 
85  T & operator() (const unsigned int i,
86  const unsigned int j);
87 
88  virtual T el(const unsigned int i,
89  const unsigned int j) const override
90  { return (*this)(i,j); }
91 
92  virtual T & el(const unsigned int i,
93  const unsigned int j) override
94  { return (*this)(i,j); }
95 
96  virtual void left_multiply (const DenseMatrixBase<T> & M2) override;
97 
101  template <typename T2>
102  void left_multiply (const DenseMatrixBase<T2> & M2);
103 
104  virtual void right_multiply (const DenseMatrixBase<T> & M2) override;
105 
109  template <typename T2>
110  void right_multiply (const DenseMatrixBase<T2> & M2);
111 
116  void vector_mult (DenseVector<T> & dest,
117  const DenseVector<T> & arg) const;
118 
124  template <typename T2>
126  const DenseVector<T2> & arg) const;
127 
133  const DenseVector<T> & arg) const;
134 
140  template <typename T2>
142  const DenseVector<T2> & arg) const;
143 
148  void vector_mult_add (DenseVector<T> & dest,
149  const T factor,
150  const DenseVector<T> & arg) const;
151 
157  template <typename T2, typename T3>
158  void vector_mult_add (DenseVector<typename CompareTypes<T, typename CompareTypes<T2,T3>::supertype>::supertype> & dest,
159  const T2 factor,
160  const DenseVector<T3> & arg) const;
161 
165  void get_principal_submatrix (unsigned int sub_m, unsigned int sub_n, DenseMatrix<T> & dest) const;
166 
170  void get_principal_submatrix (unsigned int sub_m, DenseMatrix<T> & dest) const;
171 
189  void outer_product(const DenseVector<T> & a, const DenseVector<T> & b);
190 
200  template <typename T2>
201  DenseMatrix<T> & operator = (const DenseMatrix<T2> & other_matrix);
202 
206  void swap(DenseMatrix<T> & other_matrix);
207 
212  void resize(const unsigned int new_m,
213  const unsigned int new_n);
214 
218  void scale (const T factor);
219 
223  void scale_column (const unsigned int col, const T factor);
224 
230  DenseMatrix<T> & operator *= (const T factor);
231 
237  template<typename T2, typename T3>
238  typename boostcopy::enable_if_c<
239  ScalarTraits<T2>::value, void >::type add (const T2 factor,
240  const DenseMatrix<T3> & mat);
241 
245  bool operator== (const DenseMatrix<T> & mat) const;
246 
250  bool operator!= (const DenseMatrix<T> & mat) const;
251 
258 
265 
270  Real min () const;
271 
276  Real max () const;
277 
286  Real l1_norm () const;
287 
296  Real linfty_norm () const;
297 
302 
307  template <typename T2>
309 
310 
315 
320  template <typename T2>
322 
326  T transpose (const unsigned int i,
327  const unsigned int j) const;
328 
332  void get_transpose(DenseMatrix<T> & dest) const;
333 
341  std::vector<T> & get_values() { return _val; }
342 
346  const std::vector<T> & get_values() const { return _val; }
347 
354  void condense(const unsigned int i,
355  const unsigned int j,
356  const T val,
357  DenseVector<T> & rhs)
358  { DenseMatrixBase<T>::condense (i, j, val, rhs); }
359 
365  void lu_solve (const DenseVector<T> & b,
366  DenseVector<T> & x);
367 
379  template <typename T2>
380  void cholesky_solve(const DenseVector<T2> & b,
381  DenseVector<T2> & x);
382 
391  void svd(DenseVector<Real> & sigma);
392 
404  void svd(DenseVector<Real> & sigma,
406  DenseMatrix<Number> & VT);
407 
425  void svd_solve(const DenseVector<T> & rhs,
426  DenseVector<T> & x,
427  Real rcond=std::numeric_limits<Real>::epsilon()) const;
428 
436  void evd(DenseVector<T> & lambda_real,
437  DenseVector<T> & lambda_imag);
438 
457  void evd_left(DenseVector<T> & lambda_real,
458  DenseVector<T> & lambda_imag,
459  DenseMatrix<T> & VL);
460 
479  void evd_right(DenseVector<T> & lambda_real,
480  DenseVector<T> & lambda_imag,
481  DenseMatrix<T> & VR);
482 
493  void evd_left_and_right(DenseVector<T> & lambda_real,
494  DenseVector<T> & lambda_imag,
495  DenseMatrix<T> & VL,
496  DenseMatrix<T> & VR);
497 
505  T det();
506 
518  // void inverse();
519 
526 
527 private:
528 
532  std::vector<T> _val;
533 
539  void _lu_decompose ();
540 
546  void _lu_back_substitute (const DenseVector<T> & b,
547  DenseVector<T> & x) const;
548 
556  void _cholesky_decompose();
557 
565  template <typename T2>
567  DenseVector<T2> & x) const;
568 
577 
583 
593  };
594 
602  void _multiply_blas(const DenseMatrixBase<T> & other,
603  _BLAS_Multiply_Flag flag);
604 
614  void _lu_decompose_lapack();
615 
621  void _svd_lapack(DenseVector<Real> & sigma);
622 
628  void _svd_lapack(DenseVector<Real> & sigma,
630  DenseMatrix<Number> & VT);
631 
635  void _svd_solve_lapack(const DenseVector<T> & rhs,
636  DenseVector<T> & x,
637  Real rcond) const;
638 
643  void _svd_helper (char JOBU,
644  char JOBVT,
645  std::vector<Real> & sigma_val,
646  std::vector<Number> & U_val,
647  std::vector<Number> & VT_val);
648 
657  void _evd_lapack(DenseVector<T> & lambda_real,
658  DenseVector<T> & lambda_imag,
659  DenseMatrix<T> * VL = nullptr,
660  DenseMatrix<T> * VR = nullptr);
661 
667 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
668  typedef PetscBLASInt pivot_index_t;
669 #else
670  typedef int pivot_index_t;
671 #endif
672  std::vector<pivot_index_t> _pivots;
673 
684  DenseVector<T> & x);
685 
698  void _matvec_blas(T alpha, T beta,
699  DenseVector<T> & dest,
700  const DenseVector<T> & arg,
701  bool trans=false) const;
702 };
703 
704 
705 
706 
707 
708 // ------------------------------------------------------------
712 namespace DenseMatrices
713 {
714 
720 
729 
730 }
731 
732 
733 
734 using namespace DenseMatrices;
735 
736 
737 
738 
739 
740 // ------------------------------------------------------------
741 // Dense Matrix member functions
742 template<typename T>
743 inline
744 DenseMatrix<T>::DenseMatrix(const unsigned int new_m,
745  const unsigned int new_n) :
746  DenseMatrixBase<T>(new_m,new_n),
747 #if defined(LIBMESH_HAVE_PETSC) && defined(LIBMESH_USE_REAL_NUMBERS) && defined(LIBMESH_DEFAULT_DOUBLE_PRECISION)
748  use_blas_lapack(true),
749 #else
750  use_blas_lapack(false),
751 #endif
752  _val(),
753  _decomposition_type(NONE)
754 {
755  this->resize(new_m,new_n);
756 }
757 
758 
759 
760 template<typename T>
761 inline
763 {
764  std::swap(this->_m, other_matrix._m);
765  std::swap(this->_n, other_matrix._n);
766  _val.swap(other_matrix._val);
767  DecompositionType _temp = _decomposition_type;
768  _decomposition_type = other_matrix._decomposition_type;
769  other_matrix._decomposition_type = _temp;
770 }
771 
772 
773 template <typename T>
774 template <typename T2>
775 inline
778 {
779  unsigned int mat_m = mat.m(), mat_n = mat.n();
780  this->resize(mat_m, mat_n);
781  for (unsigned int i=0; i<mat_m; i++)
782  for (unsigned int j=0; j<mat_n; j++)
783  (*this)(i,j) = mat(i,j);
784 
785  return *this;
786 }
787 
788 
789 
790 template<typename T>
791 inline
792 void DenseMatrix<T>::resize(const unsigned int new_m,
793  const unsigned int new_n)
794 {
795  _val.resize(new_m*new_n);
796 
797  this->_m = new_m;
798  this->_n = new_n;
799 
800  // zero and set decomposition_type to NONE
801  this->zero();
802 }
803 
804 
805 
806 template<typename T>
807 inline
809 {
810  _decomposition_type = NONE;
811 
812  std::fill (_val.begin(), _val.end(), static_cast<T>(0));
813 }
814 
815 
816 
817 template<typename T>
818 inline
819 T DenseMatrix<T>::operator () (const unsigned int i,
820  const unsigned int j) const
821 {
822  libmesh_assert_less (i*j, _val.size());
823  libmesh_assert_less (i, this->_m);
824  libmesh_assert_less (j, this->_n);
825 
826 
827  // return _val[(i) + (this->_m)*(j)]; // col-major
828  return _val[(i)*(this->_n) + (j)]; // row-major
829 }
830 
831 
832 
833 template<typename T>
834 inline
835 T & DenseMatrix<T>::operator () (const unsigned int i,
836  const unsigned int j)
837 {
838  libmesh_assert_less (i*j, _val.size());
839  libmesh_assert_less (i, this->_m);
840  libmesh_assert_less (j, this->_n);
841 
842  //return _val[(i) + (this->_m)*(j)]; // col-major
843  return _val[(i)*(this->_n) + (j)]; // row-major
844 }
845 
846 
847 
848 
849 
850 template<typename T>
851 inline
852 void DenseMatrix<T>::scale (const T factor)
853 {
854  for (std::size_t i=0; i<_val.size(); i++)
855  _val[i] *= factor;
856 }
857 
858 
859 template<typename T>
860 inline
861 void DenseMatrix<T>::scale_column (const unsigned int col, const T factor)
862 {
863  for (unsigned int i=0; i<this->m(); i++)
864  (*this)(i, col) *= factor;
865 }
866 
867 
868 
869 template<typename T>
870 inline
872 {
873  this->scale(factor);
874  return *this;
875 }
876 
877 
878 
879 template<typename T>
880 template<typename T2, typename T3>
881 inline
882 typename boostcopy::enable_if_c<
883  ScalarTraits<T2>::value, void >::type
884 DenseMatrix<T>::add (const T2 factor,
885  const DenseMatrix<T3> & mat)
886 {
887  libmesh_assert_equal_to (this->m(), mat.m());
888  libmesh_assert_equal_to (this->n(), mat.n());
889 
890  for (unsigned int i=0; i<this->m(); i++)
891  for (unsigned int j=0; j<this->n(); j++)
892  (*this)(i,j) += factor * mat(i,j);
893 }
894 
895 
896 
897 template<typename T>
898 inline
900 {
901  for (std::size_t i=0; i<_val.size(); i++)
902  if (_val[i] != mat._val[i])
903  return false;
904 
905  return true;
906 }
907 
908 
909 
910 template<typename T>
911 inline
913 {
914  for (std::size_t i=0; i<_val.size(); i++)
915  if (_val[i] != mat._val[i])
916  return true;
917 
918  return false;
919 }
920 
921 
922 
923 template<typename T>
924 inline
926 {
927  for (std::size_t i=0; i<_val.size(); i++)
928  _val[i] += mat._val[i];
929 
930  return *this;
931 }
932 
933 
934 
935 template<typename T>
936 inline
938 {
939  for (std::size_t i=0; i<_val.size(); i++)
940  _val[i] -= mat._val[i];
941 
942  return *this;
943 }
944 
945 
946 
947 template<typename T>
948 inline
950 {
951  libmesh_assert (this->_m);
952  libmesh_assert (this->_n);
953  Real my_min = libmesh_real((*this)(0,0));
954 
955  for (unsigned int i=0; i!=this->_m; i++)
956  {
957  for (unsigned int j=0; j!=this->_n; j++)
958  {
959  Real current = libmesh_real((*this)(i,j));
960  my_min = (my_min < current? my_min : current);
961  }
962  }
963  return my_min;
964 }
965 
966 
967 
968 template<typename T>
969 inline
971 {
972  libmesh_assert (this->_m);
973  libmesh_assert (this->_n);
974  Real my_max = libmesh_real((*this)(0,0));
975 
976  for (unsigned int i=0; i!=this->_m; i++)
977  {
978  for (unsigned int j=0; j!=this->_n; j++)
979  {
980  Real current = libmesh_real((*this)(i,j));
981  my_max = (my_max > current? my_max : current);
982  }
983  }
984  return my_max;
985 }
986 
987 
988 
989 template<typename T>
990 inline
992 {
993  libmesh_assert (this->_m);
994  libmesh_assert (this->_n);
995 
996  Real columnsum = 0.;
997  for (unsigned int i=0; i!=this->_m; i++)
998  {
999  columnsum += std::abs((*this)(i,0));
1000  }
1001  Real my_max = columnsum;
1002  for (unsigned int j=1; j!=this->_n; j++)
1003  {
1004  columnsum = 0.;
1005  for (unsigned int i=0; i!=this->_m; i++)
1006  {
1007  columnsum += std::abs((*this)(i,j));
1008  }
1009  my_max = (my_max > columnsum? my_max : columnsum);
1010  }
1011  return my_max;
1012 }
1013 
1014 
1015 
1016 template<typename T>
1017 inline
1019 {
1020  libmesh_assert (this->_m);
1021  libmesh_assert (this->_n);
1022 
1023  Real rowsum = 0.;
1024  for (unsigned int j=0; j!=this->_n; j++)
1025  {
1026  rowsum += std::abs((*this)(0,j));
1027  }
1028  Real my_max = rowsum;
1029  for (unsigned int i=1; i!=this->_m; i++)
1030  {
1031  rowsum = 0.;
1032  for (unsigned int j=0; j!=this->_n; j++)
1033  {
1034  rowsum += std::abs((*this)(i,j));
1035  }
1036  my_max = (my_max > rowsum? my_max : rowsum);
1037  }
1038  return my_max;
1039 }
1040 
1041 
1042 
1043 template<typename T>
1044 inline
1045 T DenseMatrix<T>::transpose (const unsigned int i,
1046  const unsigned int j) const
1047 {
1048  // Implement in terms of operator()
1049  return (*this)(j,i);
1050 }
1051 
1052 
1053 
1054 
1055 
1056 // template<typename T>
1057 // inline
1058 // void DenseMatrix<T>::condense(const unsigned int iv,
1059 // const unsigned int jv,
1060 // const T val,
1061 // DenseVector<T> & rhs)
1062 // {
1063 // libmesh_assert_equal_to (this->_m, rhs.size());
1064 // libmesh_assert_equal_to (iv, jv);
1065 
1066 
1067 // // move the known value into the RHS
1068 // // and zero the column
1069 // for (unsigned int i=0; i<this->m(); i++)
1070 // {
1071 // rhs(i) -= ((*this)(i,jv))*val;
1072 // (*this)(i,jv) = 0.;
1073 // }
1074 
1075 // // zero the row
1076 // for (unsigned int j=0; j<this->n(); j++)
1077 // (*this)(iv,j) = 0.;
1078 
1079 // (*this)(iv,jv) = 1.;
1080 // rhs(iv) = val;
1081 
1082 // }
1083 
1084 
1085 
1086 
1087 } // namespace libMesh
1088 
1089 
1090 
1091 
1092 #endif // LIBMESH_DENSE_MATRIX_H
T libmesh_real(T a)
double abs(double a)
void _lu_back_substitute(const DenseVector< T > &b, DenseVector< T > &x) const
virtual T & el(const unsigned int i, const unsigned int j) override
Definition: dense_matrix.h:92
T transpose(const unsigned int i, const unsigned int j) const
void swap(DenseMatrix< T > &other_matrix)
Definition: dense_matrix.h:762
void condense(const unsigned int i, const unsigned int j, const T val, DenseVector< T > &rhs)
Definition: dense_matrix.h:354
void _evd_lapack(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > *VL=nullptr, DenseMatrix< T > *VR=nullptr)
DenseMatrix< Complex > ComplexDenseMatrix
Definition: dense_matrix.h:728
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
DenseMatrix< T > & operator+=(const DenseMatrix< T > &mat)
Definition: dense_matrix.h:925
DecompositionType _decomposition_type
Definition: dense_matrix.h:582
void evd_left(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL)
if(!eq) SETERRQ2(((PetscObject) dm) -> comm, PETSC_ERR_ARG_WRONG, "DM of type %s, not of type %s",((PetscObject) dm) ->type, DMLIBMESH)
virtual ~DenseMatrix()=default
void evd_left_and_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL, DenseMatrix< T > &VR)
void _lu_back_substitute_lapack(const DenseVector< T > &b, DenseVector< T > &x)
unsigned int m() const
DenseMatrix< T > & operator-=(const DenseMatrix< T > &mat)
Definition: dense_matrix.h:937
void _svd_lapack(DenseVector< Real > &sigma)
const Number zero
Definition: libmesh.h:239
T operator()(const unsigned int i, const unsigned int j) const
Definition: dense_matrix.h:819
long double max(long double a, double b)
virtual T el(const unsigned int i, const unsigned int j) const override
Definition: dense_matrix.h:88
DenseMatrix(const unsigned int new_m=0, const unsigned int new_n=0)
Definition: dense_matrix.h:744
const std::vector< T > & get_values() const
Definition: dense_matrix.h:346
void outer_product(const DenseVector< T > &a, const DenseVector< T > &b)
void vector_mult_transpose(DenseVector< T > &dest, const DenseVector< T > &arg) const
void scale_column(const unsigned int col, const T factor)
Definition: dense_matrix.h:861
void _svd_solve_lapack(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond) const
Real l1_norm() const
Definition: dense_matrix.h:991
void _multiply_blas(const DenseMatrixBase< T > &other, _BLAS_Multiply_Flag flag)
bool operator!=(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
void _svd_helper(char JOBU, char JOBVT, std::vector< Real > &sigma_val, std::vector< Number > &U_val, std::vector< Number > &VT_val)
DenseMatrix & operator=(const DenseMatrix &)=default
void vector_mult_add(DenseVector< T > &dest, const T factor, const DenseVector< T > &arg) const
void svd(DenseVector< Real > &sigma)
void get_transpose(DenseMatrix< T > &dest) const
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseMatrix< T3 > &mat)
Definition: dense_matrix.h:884
virtual void right_multiply(const DenseMatrixBase< T > &M2) override
DenseMatrix< Real > RealDenseMatrix
Definition: dense_matrix.h:719
void evd_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VR)
void _cholesky_back_substitute(const DenseVector< T2 > &b, DenseVector< T2 > &x) const
Real linfty_norm() const
void get_principal_submatrix(unsigned int sub_m, unsigned int sub_n, DenseMatrix< T > &dest) const
std::vector< T > & get_values()
Definition: dense_matrix.h:341
bool operator==(const DenseMatrix< T > &mat) const
Definition: dense_matrix.h:899
void right_multiply_transpose(const DenseMatrix< T > &A)
void scale(const T factor)
Definition: dense_matrix.h:852
void lu_solve(const DenseVector< T > &b, DenseVector< T > &x)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void evd(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag)
void left_multiply_transpose(const DenseMatrix< T > &A)
static PetscErrorCode Mat * A
virtual void zero() override
Definition: dense_matrix.h:808
void swap(Iterator &lhs, Iterator &rhs)
void _matvec_blas(T alpha, T beta, DenseVector< T > &dest, const DenseVector< T > &arg, bool trans=false) const
void resize(const unsigned int new_m, const unsigned int new_n)
Definition: dense_matrix.h:792
bool operator==(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
Iterator & operator=(const Iterator &rhs)
virtual void left_multiply(const DenseMatrixBase< T > &M2) override
std::vector< T > _val
Definition: dense_matrix.h:532
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
DenseMatrix< T > & operator*=(const T factor)
Definition: dense_matrix.h:871
unsigned int n() const
bool operator!=(const DenseMatrix< T > &mat) const
Definition: dense_matrix.h:912
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
std::vector< pivot_index_t > _pivots
Definition: dense_matrix.h:672
long double min(long double a, double b)
void svd_solve(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond=std::numeric_limits< Real >::epsilon()) const
PetscBLASInt pivot_index_t
Definition: dense_matrix.h:668
void condense(const unsigned int i, const unsigned int j, const T val, DenseVectorBase< T > &rhs)
void vector_mult(DenseVector< T > &dest, const DenseVector< T > &arg) const