16 #ifndef dealii_full_matrix_h 17 #define dealii_full_matrix_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/numbers.h> 22 #include <deal.II/base/table.h> 23 #include <deal.II/lac/exceptions.h> 24 #include <deal.II/lac/identity_matrix.h> 25 #include <deal.II/base/tensor.h> 31 DEAL_II_NAMESPACE_OPEN
35 template <
typename number>
class Vector;
63 template <
typename number>
120 number
value()
const;
232 const number *entries);
257 template <
typename number2>
286 template <
typename number2>
296 template <
typename MatrixType>
304 template <
typename MatrixType>
317 const unsigned int src_r_i=0,
318 const unsigned int src_r_j=dim-1,
319 const unsigned int src_c_i=0,
320 const unsigned int src_c_j=dim-1,
338 const unsigned int dst_r=0,
339 const unsigned int dst_c=0)
const;
353 template <
typename MatrixType,
typename index_type>
355 const std::vector<index_type> &row_index_set,
356 const std::vector<index_type> &column_index_set);
370 template <
typename MatrixType,
typename index_type>
373 const std::vector<index_type> &column_index_set,
374 MatrixType &matrix)
const;
386 template <
typename number2>
397 template <
typename number2>
398 void fill (
const number2 *);
411 template <
typename number2>
413 const std::vector<size_type> &p_rows,
414 const std::vector<size_type> &p_cols);
480 template <
typename number2>
492 template <
typename number2>
494 const Vector<number2> &v)
const;
539 number
trace ()
const;
547 template <
class StreamType>
548 void print (StreamType &s,
549 const unsigned int width=5,
550 const unsigned int precision=2)
const;
575 const unsigned int precision=3,
576 const bool scientific =
true,
577 const unsigned int width = 0,
578 const char *zero_string =
" ",
579 const double denominator = 1.,
580 const double threshold = 0.)
const;
633 template <
typename number2>
634 void add (
const number a,
644 template <
typename number2>
645 void add (
const number a,
658 template <
typename number2>
659 void add (
const number a,
677 template <
typename number2>
690 template <
typename number2>
691 void Tadd (
const number s,
705 template <
typename number2>
729 template <
typename number2,
typename index_type>
732 const index_type *col_indices,
734 const bool elide_zero_values =
true,
735 const bool col_indices_are_sorted =
false);
788 template <
typename number2>
789 void equ (
const number a,
795 template <
typename number2>
796 void equ (
const number a,
804 template <
typename number2>
805 void equ (
const number a,
842 template <
typename number2>
853 template <
typename number2>
860 template <
typename number2>
862 const Vector<number2> &W);
869 template <
typename number2>
877 template <
typename number2>
902 template <
typename number2>
905 const bool adding=
false)
const;
925 template <
typename number2>
928 const bool adding=
false)
const;
948 template <
typename number2>
951 const bool adding=
false)
const;
972 template <
typename number2>
975 const bool adding=
false)
const;
990 const bool transpose_B =
false,
991 const bool transpose_D =
false,
992 const number scaling = number(1.));
1006 template <
typename number2>
1007 void vmult (Vector<number2> &w,
1008 const Vector<number2> &v,
1009 const bool adding=
false)
const;
1016 template <
typename number2>
1018 const Vector<number2> &v)
const;
1033 template <
typename number2>
1034 void Tvmult (Vector<number2> &w,
1035 const Vector<number2> &v,
1036 const bool adding=
false)
const;
1044 template <
typename number2>
1046 const Vector<number2> &v)
const;
1053 template <
typename somenumber>
1055 const Vector<somenumber> &src,
1056 const number omega = 1.)
const;
1064 template <
typename number2,
typename number3>
1065 number
residual (Vector<number2> &dst,
1066 const Vector<number2> &x,
1067 const Vector<number3> &b)
const;
1079 template <
typename number2>
1080 void forward (Vector<number2> &dst,
1081 const Vector<number2> &src)
const;
1090 template <
typename number2>
1091 void backward (Vector<number2> &dst,
1092 const Vector<number2> &src)
const;
1111 <<
"The maximal pivot is " << arg1
1112 <<
", which is below the threshold. The matrix may be singular.");
1118 <<
"Target region not in matrix: size in this direction=" 1119 << arg1 <<
", size of new matrix=" << arg2
1120 <<
", offset=" << arg3);
1125 "You are attempting an operation on two matrices that " 1126 "are the same object, but the operation requires that the " 1127 "two objects are in fact different.");
1145 template <
typename number>
1150 return this->n_rows();
1155 template <
typename number>
1160 return this->n_cols();
1165 template <
typename number>
1172 if (this->n_elements() != 0)
1173 this->reset_values();
1180 template <
typename number>
1181 template <
typename number2>
1190 template <
typename number>
1191 template <
typename MatrixType>
1195 this->reinit (M.m(), M.n());
1200 for (size_type row = 0; row < M.m(); ++row)
1202 const typename MatrixType::const_iterator end_row = M.end(row);
1203 for (
typename MatrixType::const_iterator entry = M.begin(row);
1204 entry != end_row; ++entry)
1205 this->el(row, entry->column()) = entry->value();
1211 template <
typename number>
1212 template <
typename MatrixType>
1216 this->reinit (M.n(), M.m());
1221 for (size_type row = 0; row < M.m(); ++row)
1223 const typename MatrixType::const_iterator end_row = M.end(row);
1224 for (
typename MatrixType::const_iterator entry = M.begin(row);
1225 entry != end_row; ++entry)
1226 this->el(entry->column(), row) = entry->value();
1232 template <
typename number>
1233 template <
typename MatrixType,
typename index_type>
1237 const std::vector<index_type> &row_index_set,
1238 const std::vector<index_type> &column_index_set)
1243 const size_type n_rows_submatrix = row_index_set.size();
1244 const size_type n_cols_submatrix = column_index_set.size();
1246 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1247 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1248 (*
this)(sub_row, sub_col) =
matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1253 template <
typename number>
1254 template <
typename MatrixType,
typename index_type>
1258 const std::vector<index_type> &column_index_set,
1259 MatrixType &matrix)
const 1264 const size_type n_rows_submatrix = row_index_set.size();
1265 const size_type n_cols_submatrix = column_index_set.size();
1267 for (size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1268 for (size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1269 matrix.set(row_index_set[sub_row],
1270 column_index_set[sub_col],
1271 (*
this)(sub_row, sub_col));
1275 template <
typename number>
1282 (*this)(i,j) = value;
1287 template <
typename number>
1288 template <
typename number2>
1291 const Vector<number2> &v)
const 1297 template <
typename number>
1298 template <
typename number2>
1301 const Vector<number2> &v)
const 1310 template <
typename number>
1323 template <
typename number>
1332 template <
typename number>
1341 template <
typename number>
1347 return matrix->el(a_row, a_col);
1351 template <
typename number>
1362 template <
typename number>
1370 if (accessor.a_col >= accessor.matrix->n())
1379 template <
typename number>
1391 template <
typename number>
1400 template <
typename number>
1409 template <
typename number>
1415 return (accessor.row() == other.accessor.row() &&
1416 accessor.column() == other.accessor.column());
1420 template <
typename number>
1426 return ! (*
this == other);
1430 template <
typename number>
1434 operator < (
const const_iterator &other)
const 1436 return (accessor.row() < other.accessor.row() ||
1437 (accessor.row() == other.accessor.row() &&
1438 accessor.column() < other.accessor.column()));
1442 template <
typename number>
1446 operator > (
const const_iterator &other)
const 1448 return (other < *
this);
1452 template <
typename number>
1457 return const_iterator(
this, 0, 0);
1461 template <
typename number>
1466 return const_iterator(
this,
m(), 0);
1470 template <
typename number>
1476 return const_iterator(
this, r, 0);
1481 template <
typename number>
1487 return const_iterator(
this, r+1, 0);
1492 template <
typename number>
1505 template <
typename number>
1506 template <
typename number2,
typename index_type>
1511 const index_type *col_indices,
1517 for (size_type col=0; col<n_cols; ++col)
1525 template <
typename number>
1526 template <
class StreamType>
1530 const unsigned int w,
1531 const unsigned int p)
const 1536 const std::streamsize old_precision = s.precision (p);
1537 const std::streamsize old_width = s.width (w);
1551 s.precision(old_precision);
1558 DEAL_II_NAMESPACE_CLOSE
types::global_dof_index size_type
number determinant() const
const Accessor & operator*() const
void diagadd(const number s)
bool operator==(const FullMatrix< number > &) const
static ::ExceptionBase & ExcEmptyMatrix()
#define AssertDimension(dim1, dim2)
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
const_iterator(const FullMatrix< number > *matrix, const size_type row, const size_type col)
number2 matrix_scalar_product(const Vector< number2 > &u, const Vector< number2 > &v) const
Contents is actually a matrix.
FullMatrix(const size_type n=0)
bool operator==(const const_iterator &) const
FullMatrix & operator/=(const number factor)
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1., const double threshold=0.) const
static ::ExceptionBase & ExcNotRegular(number arg1)
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
void cholesky(const FullMatrix< number2 > &A)
void add_row(const size_type i, const number s, const size_type j)
bool operator>(const const_iterator &) const
const_iterator begin() const
void outer_product(const Vector< number2 > &V, const Vector< number2 > &W)
#define AssertIndexRange(index, range)
void right_invert(const FullMatrix< number2 > &M)
void scatter_matrix_to(const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const
void fill_permutation(const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols)
void left_invert(const FullMatrix< number2 > &M)
void Tadd(const number s, const FullMatrix< number2 > &B)
void copy_transposed(const MatrixType &)
real_type l1_norm() const
FullMatrix & operator*=(const number factor)
const_iterator & operator++()
void set(const size_type i, const size_type j, const number value)
void triple_product(const FullMatrix< number > &A, const FullMatrix< number > &B, const FullMatrix< number > &D, const bool transpose_B=false, const bool transpose_D=false, const number scaling=number(1.))
const Accessor * operator->() const
real_type frobenius_norm() const
void copy_to(Tensor< 2, dim > &T, const size_type src_r_i=0, const size_type src_r_j=dim-1, const size_type src_c_i=0, const size_type src_c_j=dim-1, const unsigned int dst_r=0, const unsigned int dst_c=0) const
const FullMatrix< number > * matrix
void Tmmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void backward(Vector< number2 > &dst, const Vector< number2 > &src) const
numbers::NumberTraits< number >::real_type real_type
Accessor(const FullMatrix< number > *matrix, const size_type row, const size_type col)
void swap_row(const size_type i, const size_type j)
#define DeclException1(Exception1, type1, outsequence)
real_type linfty_norm() const
void invert(const FullMatrix< number2 > &M)
#define Assert(cond, exc)
bool operator!=(const const_iterator &) const
#define DeclExceptionMsg(Exception, defaulttext)
bool operator<(const const_iterator &) const
#define DeclException0(Exception0)
AlignedVector< T >::reference el(const TableIndices< N > &indices)
void extract_submatrix_from(const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set)
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
const_iterator end() const
real_type relative_symmetry_norm2() const
number residual(Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const
FullMatrix< number > & operator=(const FullMatrix< number2 > &)
void TmTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
number2 matrix_norm_square(const Vector< number2 > &v) const
void add_col(const size_type i, const number s, const size_type j)
static ::ExceptionBase & ExcIteratorPastEnd()
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
std::size_t memory_consumption() const
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void forward(Vector< number2 > &dst, const Vector< number2 > &src) const
static ::ExceptionBase & ExcSourceEqualsDestination()
void swap_col(const size_type i, const size_type j)
void copy_from(const MatrixType &)
void add(const number a, const FullMatrix< number2 > &A)
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
void print(StreamType &s, const unsigned int width=5, const unsigned int precision=2) const
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
AlignedVector< T >::size_type size_type
AlignedVector< T >::reference operator()(const TableIndices< N > &indices)
AlignedVector< T > values
void equ(const number a, const FullMatrix< number2 > &A)
#define AssertIsFinite(number)
void fill(const FullMatrix< number2 > &src, const size_type dst_offset_i=0, const size_type dst_offset_j=0, const size_type src_offset_i=0, const size_type src_offset_j=0)
void fill(InputIterator entries, const bool C_style_indexing=true)
void mTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const