|
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_chunk_sparse_matrix_h
17 # define dealii_chunk_sparse_matrix_h
36 template <
typename number>
38 template <
typename number>
53 template <
typename number,
bool Constness>
66 template <
typename number,
bool Constness>
98 template <
typename number>
148 template <
typename,
bool>
159 template <
typename number>
197 operator number()
const;
203 operator=(
const number n)
const;
209 operator+=(
const number n)
const;
215 operator-=(
const number n)
const;
221 operator*=(
const number n)
const;
227 operator/=(
const number n)
const;
279 template <
typename,
bool>
295 template <
typename number,
bool Constness>
422 template <
typename number>
681 template <
typename number2>
686 const number2 * values,
687 const bool elide_zero_values =
true,
688 const bool col_indices_are_sorted =
false);
733 template <
typename somenumber>
754 template <
typename ForwardIterator>
763 template <
typename somenumber>
778 template <
typename somenumber>
845 number * values)
const;
866 template <
class OutVector,
class InVector>
868 vmult(OutVector &dst,
const InVector &src)
const;
885 template <
class OutVector,
class InVector>
887 Tvmult(OutVector &dst,
const InVector &src)
const;
903 template <
class OutVector,
class InVector>
905 vmult_add(OutVector &dst,
const InVector &src)
const;
922 template <
class OutVector,
class InVector>
924 Tvmult_add(OutVector &dst,
const InVector &src)
const;
941 template <
typename somenumber>
948 template <
typename somenumber>
951 const Vector<somenumber> &v)
const;
959 template <
typename somenumber>
962 const Vector<somenumber> &x,
963 const Vector<somenumber> &
b)
const;
1008 template <
typename somenumber>
1011 const Vector<somenumber> &src,
1012 const number omega = 1.)
const;
1017 template <
typename somenumber>
1020 const Vector<somenumber> &src,
1021 const number om = 1.)
const;
1026 template <
typename somenumber>
1029 const Vector<somenumber> &src,
1030 const number om = 1.)
const;
1035 template <
typename somenumber>
1038 const Vector<somenumber> &src,
1039 const number om = 1.)
const;
1046 template <
typename somenumber>
1048 SSOR(Vector<somenumber> &v,
const number omega = 1.)
const;
1054 template <
typename somenumber>
1056 SOR(Vector<somenumber> &v,
const number om = 1.)
const;
1062 template <
typename somenumber>
1064 TSOR(Vector<somenumber> &v,
const number om = 1.)
const;
1076 template <
typename somenumber>
1078 PSOR(Vector<somenumber> & v,
1079 const std::vector<size_type> &permutation,
1080 const std::vector<size_type> &inverse_permutation,
1081 const number om = 1.)
const;
1093 template <
typename somenumber>
1095 TPSOR(Vector<somenumber> & v,
1096 const std::vector<size_type> &permutation,
1097 const std::vector<size_type> &inverse_permutation,
1098 const number om = 1.)
const;
1104 template <
typename somenumber>
1107 const Vector<somenumber> &
b,
1108 const number om = 1.)
const;
1114 template <
typename somenumber>
1117 const Vector<somenumber> &
b,
1118 const number om = 1.)
const;
1124 template <
typename somenumber>
1127 const Vector<somenumber> &
b,
1128 const number om = 1.)
const;
1196 begin(
const unsigned int r)
const;
1213 end(
const unsigned int r)
const;
1230 begin(
const unsigned int r);
1247 end(
const unsigned int r);
1259 print(std::ostream &out)
const;
1283 const unsigned int precision = 3,
1284 const bool scientific =
true,
1285 const unsigned int width = 0,
1286 const char * zero_string =
" ",
1287 const double denominator = 1.)
const;
1295 print_pattern(std::ostream &out,
const double threshold = 0.)
const;
1340 <<
"You are trying to access the matrix entry with index <"
1341 << arg1 <<
',' << arg2
1342 <<
">, but this entry does not exist in the sparsity pattern "
1345 "The most common cause for this problem is that you used "
1346 "a method to build the sparsity pattern that did not "
1347 "(completely) take into account all of the entries you "
1348 "will later try to write into. An example would be "
1349 "building a sparsity pattern that does not include "
1350 "the entries you will write into due to constraints "
1351 "on degrees of freedom such as hanging nodes or periodic "
1352 "boundary conditions. In such cases, building the "
1353 "sparsity pattern will succeed, but you will get errors "
1354 "such as the current one at one point or other when "
1355 "trying to write into the entries of the matrix.");
1366 <<
"The iterators denote a range of " << arg1
1367 <<
" elements, but the given number of rows was " << arg2);
1372 "You are attempting an operation on two matrices that "
1373 "are the same object, but the operation requires that the "
1374 "two objects are in fact different.");
1391 std::unique_ptr<number[]>
val;
1408 template <
typename somenumber>
1412 template <
typename,
bool>
1414 template <
typename,
bool>
1425 template <
typename number>
1434 template <
typename number>
1444 template <
typename number>
1454 template <
typename number>
1473 template <
typename number>
1484 const size_type index = compute_location(i, j);
1486 ExcInvalidIndex(i, j));
1494 template <
typename number>
1504 if (std::abs(
value) != 0.)
1506 const size_type index = compute_location(i, j);
1508 ExcInvalidIndex(i, j));
1510 val[index] +=
value;
1516 template <
typename number>
1517 template <
typename number2>
1522 const number2 * values,
1527 for (
size_type col = 0; col < n_cols; ++col)
1528 add(
row, col_indices[col],
static_cast<number
>(values[col]));
1533 template <
typename number>
1546 number * val_ptr = val.get();
1547 const number *
const end_ptr =
1550 while (val_ptr != end_ptr)
1551 *val_ptr++ *= factor;
1558 template <
typename number>
1566 const number factor_inv = 1. / factor;
1574 number * val_ptr = val.get();
1575 const number *
const end_ptr =
1579 while (val_ptr != end_ptr)
1580 *val_ptr++ *= factor_inv;
1587 template <
typename number>
1594 ExcInvalidIndex(i, j));
1595 return val[compute_location(i, j)];
1600 template <
typename number>
1605 const size_type index = compute_location(i, j);
1615 template <
typename number>
1633 template <
typename number>
1634 template <
typename ForwardIterator>
1637 const ForwardIterator
end)
1640 ExcIteratorRange(std::distance(
begin,
end), m()));
1644 using inner_iterator =
1645 typename std::iterator_traits<ForwardIterator>::value_type::const_iterator;
1647 for (ForwardIterator i =
begin; i !=
end; ++i, ++
row)
1649 const inner_iterator end_of_row = i->end();
1650 for (inner_iterator j = i->begin(); j != end_of_row; ++j)
1652 set(
row, j->first, j->second);
1663 template <
typename number>
1665 const unsigned int row)
1673 template <
typename number>
1681 template <
typename number>
1685 ,
matrix(&a.get_matrix())
1690 template <
typename number>
1695 matrix->get_sparsity_pattern().get_chunk_size();
1702 template <
typename number>
1703 inline const typename Accessor<number, true>::MatrixType &
1704 Accessor<number, true>::get_matrix()
const
1711 template <
typename number>
1712 inline Accessor<number, false>::Reference::Reference(
const Accessor *accessor,
1714 : accessor(accessor)
1718 template <
typename number>
1719 inline Accessor<number, false>::Reference::operator number()
const
1722 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1723 return accessor->matrix
1725 accessor->chunk_row *
chunk_size + accessor->chunk_col];
1730 template <
typename number>
1731 inline const typename Accessor<number, false>::Reference &
1732 Accessor<number, false>::Reference::operator=(
const number n)
const
1735 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1738 accessor->chunk_row *
chunk_size + accessor->chunk_col] = n;
1744 template <
typename number>
1745 inline const typename Accessor<number, false>::Reference &
1746 Accessor<number, false>::Reference::operator+=(
const number n)
const
1749 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1752 accessor->chunk_row *
chunk_size + accessor->chunk_col] += n;
1758 template <
typename number>
1759 inline const typename Accessor<number, false>::Reference &
1760 Accessor<number, false>::Reference::operator-=(
const number n)
const
1763 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1766 accessor->chunk_row *
chunk_size + accessor->chunk_col] -= n;
1772 template <
typename number>
1773 inline const typename Accessor<number, false>::Reference &
1774 Accessor<number, false>::Reference::operator*=(
const number n)
const
1777 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1780 accessor->chunk_row *
chunk_size + accessor->chunk_col] *= n;
1786 template <
typename number>
1787 inline const typename Accessor<number, false>::Reference &
1788 Accessor<number, false>::Reference::operator/=(
const number n)
const
1791 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1794 accessor->chunk_row *
chunk_size + accessor->chunk_col] /= n;
1800 template <
typename number>
1802 const unsigned int row)
1810 template <
typename number>
1818 template <
typename number>
1819 inline typename Accessor<number, false>::Reference
1822 return Reference(
this,
true);
1827 template <
typename number>
1828 inline typename Accessor<number, false>::MatrixType &
1829 Accessor<number, false>::get_matrix()
const
1836 template <
typename number,
bool Constness>
1838 const unsigned int row)
1844 template <
typename number,
bool Constness>
1851 template <
typename number,
bool Constness>
1859 template <
typename number,
bool Constness>
1860 inline Iterator<number, Constness> &
1868 template <
typename number,
bool Constness>
1869 inline Iterator<number, Constness>
1878 template <
typename number,
bool Constness>
1886 template <
typename number,
bool Constness>
1887 inline const Accessor<number, Constness> *Iterator<number, Constness>::
1894 template <
typename number,
bool Constness>
1898 return (accessor == other.accessor);
1902 template <
typename number,
bool Constness>
1906 return !(*
this == other);
1910 template <
typename number,
bool Constness>
1914 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1917 return (accessor < other.accessor);
1921 template <
typename number,
bool Constness>
1925 return (other < *
this);
1929 template <
typename number,
bool Constness>
1933 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1941 while (
copy != other)
1950 while (
copy != *
this)
1961 template <
typename number,
bool Constness>
1962 inline Iterator<number, Constness>
1966 for (
unsigned int i = 0; i < n; ++i)
1976 template <
typename number>
1980 return const_iterator(
this, 0);
1984 template <
typename number>
1988 return const_iterator(
this);
1992 template <
typename number>
1996 return iterator(
this, 0);
2000 template <
typename number>
2004 return iterator(
this);
2008 template <
typename number>
2013 return const_iterator(
this, r);
2018 template <
typename number>
2023 return const_iterator(
this, r + 1);
2028 template <
typename number>
2033 return iterator(
this, r);
2038 template <
typename number>
2043 return iterator(
this, r + 1);
void vmult_add(OutVector &dst, const InVector &src) const
void precondition_TSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
#define DeclExceptionMsg(Exception, defaulttext)
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
void TSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
bool operator<(const Iterator &) const
void SOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
SynchronousIterators< Iterators > operator++(SynchronousIterators< Iterators > &a)
Iterator(MatrixType *matrix, const unsigned int row)
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)
const Accessor< number, Constness > & value_type
void block_read(std::istream &in)
void SSOR(Vector< somenumber > &v, const number omega=1.) const
void block_write(std::ostream &out) const
static const bool zero_addition_can_be_elided
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
void precondition_SOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
size_type n_nonzero_elements() const
MatrixTableIterators::Accessor< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Accessor
void PSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
void print_pattern(std::ostream &out, const double threshold=0.) 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
#define AssertIndexRange(index, range)
static const size_type invalid_entry
static ::ExceptionBase & ExcDivideByZero()
const ChunkSparsityPattern & get_sparsity_pattern() const
const Accessor * accessor
Iterator operator+(const unsigned int n) const
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
number operator()(const size_type i, const size_type j) const
__global__ void set(Number *val, const Number s, const size_type N)
static ::ExceptionBase & ExcNotQuadratic()
std::unique_ptr< number[]> val
const_iterator end() const
ChunkSparseMatrix & operator*=(const number factor)
void SOR(Vector< somenumber > &v, const number om=1.) const
size_type compute_location(const size_type i, const size_type j) const
virtual ~ChunkSparseMatrix() override
Expression operator>(const Expression &lhs, const Expression &rhs)
const ChunkSparseMatrix< number > & get_matrix() const
somenumber matrix_norm_square(const Vector< somenumber > &v) const
somenumber matrix_scalar_product(const Vector< somenumber > &u, const Vector< somenumber > &v) const
ChunkSparseMatrix< number > & copy_from(const ChunkSparseMatrix< somenumber > &source)
const Accessor< number, Constness > & operator*() const
MatrixTableIterators::Iterator< TransposeTable< T >, Constness, MatrixTableIterators::Storage::column_major > Iterator
void Tvmult_add(OutVector &dst, const InVector &src) const
void Tvmult(OutVector &dst, const InVector &src) const
static const size_type invalid_entry
ChunkSparseMatrix< number > & operator=(const ChunkSparseMatrix< number > &)
#define AssertIsFinite(number)
int operator-(const Iterator &p) const
void add(const size_type i, const size_type j, const number value)
unsigned int global_dof_index
VectorType::value_type * begin(VectorType &V)
Accessor< number, Constness > accessor
#define DEAL_II_NAMESPACE_OPEN
ChunkSparseMatrix & operator/=(const number factor)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
@ matrix
Contents is actually a matrix.
size_type n_actually_nonzero_elements() const
typename numbers::NumberTraits< number >::real_type real_type
VectorType::value_type * end(VectorType &V)
real_type l1_norm() const
const Accessor< number, Constness > * operator->() const
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
static ::ExceptionBase & ExcNotInitialized()
void SSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
bool operator>(const Iterator &) const
static ::ExceptionBase & ExcSourceEqualsDestination()
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator-(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
const_iterator begin() const
types::global_dof_index size_type
real_type linfty_norm() const
void set(const size_type i, const size_type j, const number value)
typename Accessor< number, Constness >::MatrixType MatrixType
const types::global_dof_index * dummy()
static ::ExceptionBase & ExcInternalError()
types::global_dof_index size_type
friend class ChunkSparseMatrix
#define Assert(cond, exc)
void TPSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
bool operator==(const Iterator &) const
void extract_row_copy(const size_type row, const size_type array_length, size_type &row_length, size_type *column_indices, number *values) const
#define DeclException0(Exception0)
real_type frobenius_norm() const
std::size_t memory_consumption() const
number el(const size_type i, const size_type j) const
void vmult(OutVector &dst, const InVector &src) const
SmartPointer< const ChunkSparsityPattern, ChunkSparseMatrix< number > > cols
static ::ExceptionBase & ExcDifferentChunkSparsityPatterns()
number diag_element(const size_type i) const
bool operator!=(const Iterator &) const
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcIteratorRange(int arg1, int arg2)
somenumber residual(Vector< somenumber > &dst, const Vector< somenumber > &x, const Vector< somenumber > &b) const
constexpr SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator+(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertThrow(cond, exc)
virtual void reinit(const ChunkSparsityPattern &sparsity)
void precondition_SSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
void TSOR(Vector< somenumber > &v, const number om=1.) const
void copy(const T *begin, const T *end, U *dest)
Accessor(const ChunkSparsityPattern *matrix, const size_type row)
void print(std::ostream &out) const
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typename ProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)