Go to the documentation of this file.
20 #ifndef LIBMESH_DENSE_MATRIX_H
21 #define LIBMESH_DENSE_MATRIX_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/dense_matrix_base.h"
26 #include "libmesh/int_range.h"
29 #if (LIBMESH_HAVE_PETSC)
30 # include "libmesh/petsc_macro.h"
32 # define LIBMESH_SAW_I
34 # include <petscsys.h>
35 # ifndef LIBMESH_SAW_I
36 # undef I // Avoid complex.h contamination
44 #ifdef LIBMESH_HAVE_METAPHYSICL
47 template <
typename,
typename>
54 template <
typename T,
typename D>
56 template <
typename T,
typename D>
65 template <
typename T>
class DenseVector;
78 class DenseMatrix :
public DenseMatrixBase<T>
86 const unsigned int new_n=0);
98 virtual void zero()
override;
104 unsigned int col_id,
unsigned int col_size)
const;
110 const unsigned int j)
const;
116 const unsigned int j);
118 virtual T
el(
const unsigned int i,
119 const unsigned int j)
const override
120 {
return (*
this)(i,j); }
122 virtual T &
el(
const unsigned int i,
123 const unsigned int j)
override
124 {
return (*
this)(i,j); }
131 template <
typename T2>
139 template <
typename T2>
154 template <
typename T2>
170 template <
typename T2>
187 template <
typename T2,
typename T3>
230 template <
typename T2>
242 void resize(
const unsigned int new_m,
243 const unsigned int new_n);
248 void scale (
const T factor);
253 void scale_column (
const unsigned int col,
const T factor);
267 template<
typename T2,
typename T3>
337 template <typename T2>
350 template <typename T2>
357 const
unsigned int j) const;
385 const unsigned int j,
409 template <
typename T2>
457 Real rcond=std::numeric_limits<Real>::epsilon())
const;
595 template <
typename T2>
675 std::vector<Real> & sigma_val,
676 std::vector<Number> & U_val,
677 std::vector<Number> & VT_val);
697 #if (LIBMESH_HAVE_PETSC && LIBMESH_USE_REAL_NUMBERS)
731 bool trans=
false)
const;
742 namespace DenseMatrices
764 using namespace DenseMatrices;
775 const unsigned int new_n) :
777 #
if defined(LIBMESH_HAVE_PETSC) && defined(LIBMESH_USE_REAL_NUMBERS) && defined(LIBMESH_DEFAULT_DOUBLE_PRECISION)
778 use_blas_lapack(true),
780 use_blas_lapack(false),
783 _decomposition_type(NONE)
785 this->resize(new_m,new_n);
796 _val.swap(other_matrix._val);
797 DecompositionType _temp = _decomposition_type;
798 _decomposition_type = other_matrix._decomposition_type;
799 other_matrix._decomposition_type = _temp;
803 template <
typename T>
804 template <
typename T2>
809 unsigned int mat_m = mat.m(), mat_n = mat.n();
810 this->resize(mat_m, mat_n);
811 for (
unsigned int i=0; i<mat_m; i++)
812 for (
unsigned int j=0; j<mat_n; j++)
813 (*
this)(i,j) = mat(i,j);
823 const unsigned int new_n)
840 _decomposition_type = NONE;
842 std::fill (_val.begin(), _val.end(), static_cast<T>(0));
850 unsigned int col_id,
unsigned int col_size)
const
852 libmesh_assert_less (row_id + row_size - 1, this->_m);
853 libmesh_assert_less (col_id + col_size - 1, this->_n);
858 sub._val.resize(row_size * col_size);
860 unsigned int end_col = this->_n - col_size - col_id;
861 unsigned int p = row_id * this->_n;
863 for (
unsigned int i=0; i<row_size; i++)
867 for (
unsigned int j=0; j<col_size; j++)
868 sub._val[q++] = _val[p++];
881 const unsigned int j)
const
883 libmesh_assert_less (i*j, _val.size());
884 libmesh_assert_less (i, this->_m);
885 libmesh_assert_less (j, this->_n);
889 return _val[(i)*(this->_n) + (j)];
897 const unsigned int j)
899 libmesh_assert_less (i*j, _val.size());
900 libmesh_assert_less (i, this->_m);
901 libmesh_assert_less (j, this->_n);
904 return _val[(i)*(this->_n) + (j)];
915 for (
auto & v : _val)
925 (*this)(i, col) *= factor;
941 template<
typename T2,
typename T3>
948 libmesh_assert_equal_to (this->m(), mat.m());
949 libmesh_assert_equal_to (this->n(), mat.n());
953 (*
this)(i,j) += factor * mat(i,j);
963 if (_val[i] != mat._val[i])
976 if (_val[i] != mat._val[i])
989 _val[i] += mat._val[i];
1001 _val[i] -= mat._val[i];
1008 template<
typename T>
1016 for (
unsigned int i=0; i!=this->_m; i++)
1018 for (
unsigned int j=0; j!=this->_n; j++)
1021 my_min = (my_min < current? my_min : current);
1029 template<
typename T>
1037 for (
unsigned int i=0; i!=this->_m; i++)
1039 for (
unsigned int j=0; j!=this->_n; j++)
1042 my_max = (my_max > current? my_max : current);
1050 template<
typename T>
1058 for (
unsigned int i=0; i!=this->_m; i++)
1060 columnsum +=
std::abs((*
this)(i,0));
1062 auto my_max = columnsum;
1063 for (
unsigned int j=1; j!=this->_n; j++)
1066 for (
unsigned int i=0; i!=this->_m; i++)
1068 columnsum +=
std::abs((*
this)(i,j));
1070 my_max = (my_max > columnsum? my_max : columnsum);
1077 template<
typename T>
1085 for (
unsigned int j=0; j!=this->_n; j++)
1089 auto my_max = rowsum;
1090 for (
unsigned int i=1; i!=this->_m; i++)
1093 for (
unsigned int j=0; j!=this->_n; j++)
1097 my_max = (my_max > rowsum? my_max : rowsum);
1104 template<
typename T>
1107 const unsigned int j)
const
1110 return (*
this)(j,i);
1153 #endif // LIBMESH_DENSE_MATRIX_H
void _svd_lapack(DenseVector< Real > &sigma)
Computes an SVD of the matrix using the Lapack routine "getsvd".
auto min() const -> decltype(libmesh_real(T(0)))
bool use_blas_lapack
Computes the inverse of the dense matrix (assuming it is invertible) by first computing the LU decomp...
void get_transpose(DenseMatrix< T > &dest) const
Put the tranposed matrix into dest.
DenseMatrix< Complex > ComplexDenseMatrix
This typedef may be either a real-only matrix, or a truly complex matrix, depending on how Number was...
std::vector< pivot_index_t > _pivots
void left_multiply_transpose(const DenseMatrix< T > &A)
Left multiplies by the transpose of the matrix A.
const std::vector< T > & get_values() const
void right_multiply_transpose(const DenseMatrix< T > &A)
Right multiplies by the transpose of the matrix A.
bool operator==(const DenseMatrix< T > &mat) const
void _svd_solve_lapack(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond) const
Called by svd_solve(rhs).
void _cholesky_back_substitute(const DenseVector< T2 > &b, DenseVector< T2 > &x) const
Solves the equation Ax=b for the unknown value x and rhs b based on the Cholesky factorization of A.
std::vector< T > _val
The actual data values, stored as a 1D array.
DenseMatrix(const unsigned int new_m=0, const unsigned int new_n=0)
Constructor.
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
The libMesh namespace provides an interface to certain functionality in the library.
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add(const T2 factor, const DenseMatrix< T3 > &mat)
Adds factor times mat to this matrix.
void _svd_helper(char JOBU, char JOBVT, std::vector< Real > &sigma_val, std::vector< Number > &U_val, std::vector< Number > &VT_val)
Helper function that actually performs the SVD.
void evd_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VR)
Compute the eigenvalues (both real and imaginary parts) and right eigenvectors of a general matrix,...
void evd(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag)
Compute the eigenvalues (both real and imaginary parts) of a general matrix.
Defines a dense matrix for use in Finite Element-type computations.
void _cholesky_decompose()
Decomposes a symmetric positive definite matrix into a product of two lower triangular matrices accor...
void scale_column(const unsigned int col, const T factor)
Multiplies every element in the column col matrix by factor.
DenseMatrix< T > & operator*=(const T factor)
Multiplies every element in the matrix by factor.
std::vector< T > & get_values()
DenseMatrix sub_matrix(unsigned int row_id, unsigned int row_size, unsigned int col_id, unsigned int col_size) const
Get submatrix with the smallest row and column indices and the submatrix size.
void condense(const unsigned int i, const unsigned int j, const T val, DenseVectorBase< T > &rhs)
Condense-out the (i,j) entry of the matrix, forcing it to take on the value val.
virtual void right_multiply(const DenseMatrixBase< T > &M2) override
Performs the operation: (*this) <- (*this) * M3.
void resize(const unsigned int new_m, const unsigned int new_n)
Resize the matrix.
DenseMatrix< T > & operator+=(const DenseMatrix< T > &mat)
Adds mat to this matrix.
auto max() const -> decltype(libmesh_real(T(0)))
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
void swap(DenseMatrix< T > &other_matrix)
STL-like swap method.
void vector_mult(DenseVector< T > &dest, const DenseVector< T > &arg) const
Performs the matrix-vector multiplication, dest := (*this) * arg.
virtual ~DenseMatrix()=default
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
void _matvec_blas(T alpha, T beta, DenseVector< T > &dest, const DenseVector< T > &arg, bool trans=false) const
Uses the BLAS GEMV function (through PETSc) to compute.
Defines an abstract dense matrix base class for use in Finite Element-type computations.
bool operator!=(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
PetscBLASInt pivot_index_t
Array used to store pivot indices.
void _lu_back_substitute(const DenseVector< T > &b, DenseVector< T > &x) const
Solves the system Ax=b through back substitution.
void vector_mult_transpose(DenseVector< T > &dest, const DenseVector< T > &arg) const
Performs the matrix-vector multiplication, dest := (*this)^T * arg.
MetaPhysicL::DualNumber< T, D > abs(MetaPhysicL::DualNumber< T, D > &&in)
void evd_left_and_right(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL, DenseMatrix< T > &VR)
Compute the eigenvalues (both real and imaginary parts) as well as the left and right eigenvectors of...
void vector_mult_add(DenseVector< T > &dest, const T factor, const DenseVector< T > &arg) const
Performs the scaled matrix-vector multiplication, dest += factor * (*this) * arg.
static PetscErrorCode Mat * A
void svd_solve(const DenseVector< T > &rhs, DenseVector< T > &x, Real rcond=std::numeric_limits< Real >::epsilon()) const
Solve the system of equations for in the least-squares sense.
_BLAS_Multiply_Flag
Enumeration used to determine the behavior of the _multiply_blas function.
void condense(const unsigned int i, const unsigned int j, const T val, DenseVector< T > &rhs)
Condense-out the (i,j) entry of the matrix, forcing it to take on the value val.
void evd_left(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > &VL)
Compute the eigenvalues (both real and imaginary parts) and left eigenvectors of a general matrix,...
virtual void zero() override
Set every element in the matrix to 0.
virtual T & el(const unsigned int i, const unsigned int j) override
void _lu_back_substitute_lapack(const DenseVector< T > &b, DenseVector< T > &x)
Companion function to _lu_decompose_lapack().
void cholesky_solve(const DenseVector< T2 > &b, DenseVector< T2 > &x)
For symmetric positive definite (SPD) matrices.
void _lu_decompose_lapack()
Computes an LU factorization of the matrix using the Lapack routine "getrf".
bool operator!=(const DenseMatrix< T > &mat) const
DenseMatrix< Real > RealDenseMatrix
Convenient definition of a real-only dense matrix.
void swap(Iterator &lhs, Iterator &rhs)
swap, used to implement op=
virtual void left_multiply(const DenseMatrixBase< T > &M2) override
Performs the operation: (*this) <- M2 * (*this)
DecompositionType
The decomposition schemes above change the entries of the matrix A.
void _multiply_blas(const DenseMatrixBase< T > &other, _BLAS_Multiply_Flag flag)
The _multiply_blas function computes A <- op(A) * op(B) using BLAS gemm function.
Iterator & operator=(const Iterator &rhs)
Assignment operator.
void get_principal_submatrix(unsigned int sub_m, unsigned int sub_n, DenseMatrix< T > &dest) const
Put the sub_m x sub_n principal submatrix into dest.
void outer_product(const DenseVector< T > &a, const DenseVector< T > &b)
Computes the outer (dyadic) product of two vectors and stores in (*this).
void svd(DenseVector< Real > &sigma)
Compute the singular value decomposition of the matrix.
T transpose(const unsigned int i, const unsigned int j) const
T operator()(const unsigned int i, const unsigned int j) const
void _evd_lapack(DenseVector< T > &lambda_real, DenseVector< T > &lambda_imag, DenseMatrix< T > *VL=nullptr, DenseMatrix< T > *VR=nullptr)
Computes the eigenvalues of the matrix using the Lapack routine "DGEEV".
bool operator==(const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &x, const variant_filter_iterator< Predicate, Type, ReferenceType, PointerType > &y)
auto l1_norm() const -> decltype(std::abs(T(0)))
virtual T el(const unsigned int i, const unsigned int j) const override
void lu_solve(const DenseVector< T > &b, DenseVector< T > &x)
Solve the system Ax=b given the input vector b.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void _lu_decompose()
Form the LU decomposition of the matrix.
void scale(const T factor)
Multiplies every element in the matrix by factor.
DenseMatrix & operator=(const DenseMatrix &)=default
Defines a dense vector for use in Finite Element-type computations.
auto linfty_norm() const -> decltype(std::abs(T(0)))
DecompositionType _decomposition_type
This flag keeps track of which type of decomposition has been performed on the matrix.
DenseMatrix< T > & operator-=(const DenseMatrix< T > &mat)
Subtracts mat from this matrix.