|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
16 #ifndef dealii_full_matrix_h
17 #define dealii_full_matrix_h
40 template <
typename number>
42 template <
typename number>
70 template <
typename number>
78 "The FullMatrix class does not support auto-differentiable numbers.");
172 template <
typename number2>
201 template <
typename number2>
211 template <
typename MatrixType>
220 template <
typename MatrixType>
234 const unsigned int src_r_i = 0,
235 const unsigned int src_r_j = dim - 1,
236 const unsigned int src_c_i = 0,
237 const unsigned int src_c_j = dim - 1,
254 const unsigned int dst_r = 0,
255 const unsigned int dst_c = 0)
const;
269 template <
typename MatrixType,
typename index_type>
272 const std::vector<index_type> &row_index_set,
273 const std::vector<index_type> &column_index_set);
287 template <
typename MatrixType,
typename index_type>
290 const std::vector<index_type> &column_index_set,
291 MatrixType &
matrix)
const;
303 template <
typename number2>
315 template <
typename number2>
317 fill(
const number2 *);
330 template <
typename number2>
333 const std::vector<size_type> &p_rows,
334 const std::vector<size_type> &p_cols);
403 template <
typename number2>
416 template <
typename number2>
419 const Vector<number2> &v)
const;
478 template <
class StreamType>
480 print(StreamType & s,
481 const unsigned int width = 5,
482 const unsigned int precision = 2)
const;
508 const unsigned int precision = 3,
509 const bool scientific =
true,
510 const unsigned int width = 0,
511 const char * zero_string =
" ",
512 const double denominator = 1.,
513 const double threshold = 0.)
const;
573 template <
typename number2>
584 template <
typename number2>
599 template <
typename number2>
619 template <
typename number2>
633 template <
typename number2>
648 template <
typename number2>
672 template <
typename number2,
typename index_type>
676 const index_type *col_indices,
678 const bool elide_zero_values =
true,
679 const bool col_indices_are_sorted =
false);
737 template <
typename number2>
744 template <
typename number2>
754 template <
typename number2>
795 template <
typename number2>
807 template <
typename number2>
815 template <
typename number2>
817 outer_product(
const Vector<number2> &
V,
const Vector<number2> &W);
824 template <
typename number2>
833 template <
typename number2>
859 template <
typename number2>
863 const bool adding =
false)
const;
883 template <
typename number2>
887 const bool adding =
false)
const;
907 template <
typename number2>
911 const bool adding =
false)
const;
932 template <
typename number2>
936 const bool adding =
false)
const;
952 const bool transpose_B =
false,
953 const bool transpose_D =
false,
954 const number scaling = number(1.));
968 template <
typename number2>
970 vmult(Vector<number2> &
w,
971 const Vector<number2> &v,
972 const bool adding =
false)
const;
979 template <
typename number2>
981 vmult_add(Vector<number2> &
w,
const Vector<number2> &v)
const;
996 template <
typename number2>
999 const Vector<number2> &v,
1000 const bool adding =
false)
const;
1008 template <
typename number2>
1010 Tvmult_add(Vector<number2> &
w,
const Vector<number2> &v)
const;
1017 template <
typename somenumber>
1020 const Vector<somenumber> &src,
1021 const number omega = 1.)
const;
1029 template <
typename number2,
typename number3>
1032 const Vector<number2> &x,
1033 const Vector<number3> &
b)
const;
1045 template <
typename number2>
1047 forward(Vector<number2> &dst,
const Vector<number2> &src)
const;
1056 template <
typename number2>
1058 backward(Vector<number2> &dst,
const Vector<number2> &src)
const;
1078 <<
"The maximal pivot is " << arg1
1079 <<
", which is below the threshold. The matrix may be singular.");
1087 <<
"Target region not in matrix: size in this direction="
1088 << arg1 <<
", size of new matrix=" << arg2
1089 <<
", offset=" << arg3);
1094 "You are attempting an operation on two matrices that "
1095 "are the same object, but the operation requires that the "
1096 "two objects are in fact different.");
1111 template <
typename number>
1115 return this->n_rows();
1120 template <
typename number>
1124 return this->n_cols();
1129 template <
typename number>
1136 if (this->n_elements() != 0)
1137 this->reset_values();
1144 template <
typename number>
1145 template <
typename number2>
1154 template <
typename number>
1155 template <
typename MatrixType>
1159 this->
reinit(M.m(), M.n());
1164 for (
size_type row = 0; row < M.m(); ++row)
1166 const typename MatrixType::const_iterator end_row = M.end(row);
1167 for (
typename MatrixType::const_iterator entry = M.begin(row);
1170 this->el(row, entry->column()) = entry->value();
1176 template <
typename number>
1177 template <
typename MatrixType>
1181 this->
reinit(M.n(), M.m());
1186 for (
size_type row = 0; row < M.m(); ++row)
1188 const typename MatrixType::const_iterator end_row = M.end(row);
1189 for (
typename MatrixType::const_iterator entry = M.begin(row);
1192 this->el(entry->column(), row) = entry->value();
1198 template <
typename number>
1199 template <
typename MatrixType,
typename index_type>
1202 const MatrixType &
matrix,
1203 const std::vector<index_type> &row_index_set,
1204 const std::vector<index_type> &column_index_set)
1209 const size_type n_rows_submatrix = row_index_set.size();
1210 const size_type n_cols_submatrix = column_index_set.size();
1212 for (
size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1213 for (
size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1214 (*
this)(sub_row, sub_col) =
1215 matrix.el(row_index_set[sub_row], column_index_set[sub_col]);
1220 template <
typename number>
1221 template <
typename MatrixType,
typename index_type>
1224 const std::vector<index_type> &row_index_set,
1225 const std::vector<index_type> &column_index_set,
1226 MatrixType &
matrix)
const
1231 const size_type n_rows_submatrix = row_index_set.size();
1232 const size_type n_cols_submatrix = column_index_set.size();
1234 for (
size_type sub_row = 0; sub_row < n_rows_submatrix; ++sub_row)
1235 for (
size_type sub_col = 0; sub_col < n_cols_submatrix; ++sub_col)
1236 matrix.set(row_index_set[sub_row],
1237 column_index_set[sub_col],
1238 (*
this)(sub_row, sub_col));
1242 template <
typename number>
1248 (*this)(i, j) =
value;
1253 template <
typename number>
1254 template <
typename number2>
1257 const Vector<number2> &v)
const
1263 template <
typename number>
1264 template <
typename number2>
1267 const Vector<number2> &v)
const
1274 template <
typename number>
1279 return iterator(
this, r, 0);
1284 template <
typename number>
1289 return iterator(
this, r + 1, 0);
1294 template <
typename number>
1299 return const_iterator(
this, r, 0);
1304 template <
typename number>
1309 return const_iterator(
this, r + 1, 0);
1314 template <
typename number>
1321 this->operator()(r, c) += v;
1326 template <
typename number>
1327 template <
typename number2,
typename index_type>
1331 const index_type *col_indices,
1332 const number2 * values,
1337 for (
size_type col = 0; col < n_cols; ++col)
1340 this->operator()(row, col_indices[col]) += values[col];
1345 template <
typename number>
1346 template <
class StreamType>
1349 const unsigned int w,
1350 const unsigned int p)
const
1352 Assert(!this->empty(), ExcEmptyMatrix());
1355 const std::streamsize old_precision = s.precision(p);
1356 const std::streamsize old_width = s.width(
w);
1358 for (
size_type i = 0; i < this->m(); ++i)
1360 for (
size_type j = 0; j < this->n(); ++j)
1364 s << this->el(i, j);
1370 s.precision(old_precision);
iterator end(const size_type r)
void Tvmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
#define DeclExceptionMsg(Exception, defaulttext)
void TmTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void equ(const number a, const FullMatrix< number2 > &A)
FullMatrix< number > & operator=(const FullMatrix< number2 > &)
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
void copy_from(const MatrixType &)
static ::ExceptionBase & ExcNotRegular(number arg1)
static ::ExceptionBase & ExcEmptyMatrix()
void copy_transposed(const MatrixType &)
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
void forward(Vector< number2 > &dst, const Vector< number2 > &src) const
void add(const number a, const FullMatrix< number2 > &A)
FullMatrix & operator*=(const number factor)
void Tmmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void set(const size_type i, const size_type j, const number value)
number determinant() const
void add_col(const size_type i, const number s, const size_type j)
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)
static ::ExceptionBase & ExcSourceEqualsDestination()
real_type relative_symmetry_norm2() const
#define AssertIndexRange(index, range)
SymmetricTensor< 2, dim, Number > C(const Tensor< 2, dim, Number > &F)
real_type linfty_norm() const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void print(StreamType &s, const unsigned int width=5, const unsigned int precision=2) const
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
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
iterator begin(const size_type r)
number residual(Vector< number2 > &dst, const Vector< number2 > &x, const Vector< number3 > &b) const
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
void mTmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
typename AlignedVector< CoefficientType >::size_type size_type
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
#define DEAL_II_NAMESPACE_OPEN
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
@ matrix
Contents is actually a matrix.
typename Table< 2, CoefficientType >::iterator iterator
void swap_col(const size_type i, const size_type j)
#define DeclException3(Exception3, type1, type2, type3, outsequence)
number2 matrix_scalar_product(const Vector< number2 > &u, const Vector< number2 > &v) const
void fill(InputIterator entries, const bool C_style_indexing=true)
void outer_product(const Vector< number2 > &V, const Vector< number2 > &W)
FullMatrix(const size_type n=0)
#define DeclException1(Exception1, type1, outsequence)
#define AssertDimension(dim1, dim2)
void diagadd(const number s)
void cholesky(const FullMatrix< number2 > &A)
void swap_row(const size_type i, const size_type j)
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.))
typename numbers::NumberTraits< CoefficientType >::real_type real_type
static ::ExceptionBase & ExcMatrixNotPositiveDefinite()
real_type frobenius_norm() const
static ::ExceptionBase & ExcInvalidDestination(size_type arg1, size_type arg2, size_type arg3)
types::global_dof_index size_type
#define Assert(cond, exc)
void Tadd(const number s, const FullMatrix< number2 > &B)
#define DeclException0(Exception0)
void backward(Vector< number2 > &dst, const Vector< number2 > &src) const
void scatter_matrix_to(const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set, MatrixType &matrix) const
void mmult(FullMatrix< number2 > &C, const FullMatrix< number2 > &B, const bool adding=false) const
void invert(const FullMatrix< number2 > &M)
void left_invert(const FullMatrix< number2 > &M)
typename Table< 2, CoefficientType >::const_iterator const_iterator
static ::ExceptionBase & ExcScalarAssignmentOnlyForZeroValue()
std::size_t memory_consumption() const
void right_invert(const FullMatrix< number2 > &M)
#define DEAL_II_NAMESPACE_CLOSE
void Tvmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
bool operator==(const FullMatrix< number > &) const
AlignedVector< number > values
CoefficientType value_type
void add_row(const size_type i, const number s, const size_type j)
void vmult_add(Vector< number2 > &w, const Vector< number2 > &v) const
real_type l1_norm() const
FullMatrix & operator/=(const number factor)
void fill_permutation(const FullMatrix< number2 > &src, const std::vector< size_type > &p_rows, const std::vector< size_type > &p_cols)
void extract_submatrix_from(const MatrixType &matrix, const std::vector< index_type > &row_index_set, const std::vector< index_type > &column_index_set)
number2 matrix_norm_square(const Vector< number2 > &v) const