16#ifndef dealii_chunk_sparse_matrix_h
17#define dealii_chunk_sparse_matrix_h
37template <
typename number>
39template <
typename number>
55 template <
typename number,
bool Constness>
68 template <
typename number,
bool Constness>
100 template <
typename number>
150 template <
typename,
bool>
161 template <
typename number>
199 operator number()
const;
281 template <
typename,
bool>
297 template <
typename number,
bool Constness>
416 template <
typename number,
bool Constness>
422 Iterator<number, Constness>::value_type;
424 Iterator<number, Constness>::difference_type;
448template <
typename number>
707 template <
typename number2>
712 const number2 * values,
713 const bool elide_zero_values =
true,
714 const bool col_indices_are_sorted =
false);
759 template <
typename somenumber>
780 template <
typename ForwardIterator>
789 template <
typename somenumber>
804 template <
typename somenumber>
871 number * values)
const;
892 template <
class OutVector,
class InVector>
894 vmult(OutVector &dst,
const InVector &src)
const;
911 template <
class OutVector,
class InVector>
913 Tvmult(OutVector &dst,
const InVector &src)
const;
929 template <
class OutVector,
class InVector>
948 template <
class OutVector,
class InVector>
967 template <
typename somenumber>
974 template <
typename somenumber>
985 template <
typename somenumber>
1034 template <
typename somenumber>
1038 const number omega = 1.)
const;
1043 template <
typename somenumber>
1047 const number om = 1.)
const;
1052 template <
typename somenumber>
1056 const number om = 1.)
const;
1061 template <
typename somenumber>
1065 const number om = 1.)
const;
1072 template <
typename somenumber>
1080 template <
typename somenumber>
1088 template <
typename somenumber>
1102 template <
typename somenumber>
1105 const std::vector<size_type> &permutation,
1106 const std::vector<size_type> &inverse_permutation,
1107 const number om = 1.)
const;
1119 template <
typename somenumber>
1122 const std::vector<size_type> &permutation,
1123 const std::vector<size_type> &inverse_permutation,
1124 const number om = 1.)
const;
1130 template <
typename somenumber>
1134 const number om = 1.)
const;
1140 template <
typename somenumber>
1144 const number om = 1.)
const;
1150 template <
typename somenumber>
1154 const number om = 1.)
const;
1239 end(
const unsigned int r)
const;
1309 const unsigned int precision = 3,
1310 const bool scientific =
true,
1311 const unsigned int width = 0,
1312 const char * zero_string =
" ",
1313 const double denominator = 1.)
const;
1366 <<
"You are trying to access the matrix entry with index <"
1367 << arg1 <<
',' << arg2
1368 <<
">, but this entry does not exist in the sparsity pattern "
1371 "The most common cause for this problem is that you used "
1372 "a method to build the sparsity pattern that did not "
1373 "(completely) take into account all of the entries you "
1374 "will later try to write into. An example would be "
1375 "building a sparsity pattern that does not include "
1376 "the entries you will write into due to constraints "
1377 "on degrees of freedom such as hanging nodes or periodic "
1378 "boundary conditions. In such cases, building the "
1379 "sparsity pattern will succeed, but you will get errors "
1380 "such as the current one at one point or other when "
1381 "trying to write into the entries of the matrix.");
1392 <<
"The iterators denote a range of " << arg1
1393 <<
" elements, but the given number of rows was " << arg2);
1398 "You are attempting an operation on two vectors that "
1399 "are the same object, but the operation requires that the "
1400 "two objects are in fact different.");
1417 std::unique_ptr<number[]>
val;
1434 template <
typename somenumber>
1438 template <
typename,
bool>
1440 template <
typename,
bool>
1451template <
typename number>
1460template <
typename number>
1470template <
typename number>
1480template <
typename number>
1487 cols->sparsity_pattern(i / chunk_size, j / chunk_size);
1493 return (chunk_index * chunk_size * chunk_size +
1494 (i % chunk_size) * chunk_size + (j % chunk_size));
1499template <
typename number>
1512 ExcInvalidIndex(i, j));
1520template <
typename number>
1534 ExcInvalidIndex(i, j));
1542template <
typename number>
1543template <
typename number2>
1548 const number2 * values,
1553 for (
size_type col = 0; col < n_cols; ++col)
1554 add(
row, col_indices[col],
static_cast<number
>(values[col]));
1559template <
typename number>
1572 number * val_ptr = val.get();
1573 const number *
const end_ptr =
1576 while (val_ptr != end_ptr)
1577 *val_ptr++ *= factor;
1584template <
typename number>
1592 const number factor_inv = 1. / factor;
1600 number * val_ptr = val.get();
1601 const number *
const end_ptr =
1605 while (val_ptr != end_ptr)
1606 *val_ptr++ *= factor_inv;
1613template <
typename number>
1620 ExcInvalidIndex(i, j));
1621 return val[compute_location(i, j)];
1626template <
typename number>
1641template <
typename number>
1654 (i %
chunk_size) * chunk_size + (i % chunk_size)];
1659template <
typename number>
1660template <
typename ForwardIterator>
1663 const ForwardIterator end)
1666 ExcIteratorRange(std::distance(begin, end), m()));
1670 using inner_iterator =
1671 typename std::iterator_traits<ForwardIterator>::value_type::const_iterator;
1673 for (ForwardIterator i = begin; i !=
end; ++i, ++
row)
1675 const inner_iterator end_of_row = i->end();
1676 for (inner_iterator j = i->begin(); j != end_of_row; ++j)
1678 set(
row, j->first, j->second);
1689 template <
typename number>
1691 const unsigned int row)
1699 template <
typename number>
1700 inline Accessor<number, true>::Accessor(
const MatrixType *matrix)
1707 template <
typename number>
1708 inline Accessor<number, true>::Accessor(
1711 ,
matrix(&a.get_matrix())
1716 template <
typename number>
1718 Accessor<number, true>::value()
const
1721 matrix->get_sparsity_pattern().get_chunk_size();
1728 template <
typename number>
1729 inline const typename Accessor<number, true>::MatrixType &
1730 Accessor<number, true>::get_matrix()
const
1737 template <
typename number>
1738 inline Accessor<number, false>::Reference::Reference(
const Accessor *accessor,
1740 : accessor(accessor)
1744 template <
typename number>
1745 inline Accessor<number, false>::Reference::operator number()
const
1748 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1749 return accessor->matrix
1751 accessor->chunk_row *
chunk_size + accessor->chunk_col];
1756 template <
typename number>
1757 inline const typename Accessor<number, false>::Reference &
1758 Accessor<number, false>::Reference::operator=(
const number n)
const
1761 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1764 accessor->chunk_row *
chunk_size + accessor->chunk_col] = n;
1770 template <
typename number>
1771 inline const typename Accessor<number, false>::Reference &
1772 Accessor<number, false>::Reference::operator+=(
const number n)
const
1775 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1778 accessor->chunk_row *
chunk_size + accessor->chunk_col] += n;
1784 template <
typename number>
1785 inline const typename Accessor<number, false>::Reference &
1786 Accessor<number, false>::Reference::operator-=(
const number n)
const
1789 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1792 accessor->chunk_row *
chunk_size + accessor->chunk_col] -= n;
1798 template <
typename number>
1799 inline const typename Accessor<number, false>::Reference &
1800 Accessor<number, false>::Reference::operator*=(
const number n)
const
1803 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1806 accessor->chunk_row *
chunk_size + accessor->chunk_col] *= n;
1812 template <
typename number>
1813 inline const typename Accessor<number, false>::Reference &
1814 Accessor<number, false>::Reference::operator/=(
const number n)
const
1817 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1820 accessor->chunk_row *
chunk_size + accessor->chunk_col] /= n;
1826 template <
typename number>
1827 inline Accessor<number, false>::Accessor(MatrixType * matrix,
1828 const unsigned int row)
1836 template <
typename number>
1837 inline Accessor<number, false>::Accessor(MatrixType *matrix)
1844 template <
typename number>
1845 inline typename Accessor<number, false>::Reference
1846 Accessor<number, false>::value()
const
1848 return Reference(
this,
true);
1853 template <
typename number>
1854 inline typename Accessor<number, false>::MatrixType &
1855 Accessor<number, false>::get_matrix()
const
1862 template <
typename number,
bool Constness>
1863 inline Iterator<number, Constness>::Iterator(MatrixType * matrix,
1864 const unsigned int row)
1870 template <
typename number,
bool Constness>
1871 inline Iterator<number, Constness>::Iterator(MatrixType *matrix)
1877 template <
typename number,
bool Constness>
1878 inline Iterator<number, Constness>::Iterator(
1885 template <
typename number,
bool Constness>
1886 inline Iterator<number, Constness> &
1887 Iterator<number, Constness>::operator++()
1894 template <
typename number,
bool Constness>
1895 inline Iterator<number, Constness>
1896 Iterator<number, Constness>::operator++(
int)
1898 const Iterator iter = *
this;
1904 template <
typename number,
bool Constness>
1905 inline const Accessor<number, Constness> &
1906 Iterator<number, Constness>::operator*()
const
1912 template <
typename number,
bool Constness>
1913 inline const Accessor<number, Constness> *
1914 Iterator<number, Constness>::operator->()
const
1920 template <
typename number,
bool Constness>
1922 Iterator<number, Constness>::operator==(
const Iterator &other)
const
1924 return (accessor == other.accessor);
1928 template <
typename number,
bool Constness>
1930 Iterator<number, Constness>::operator!=(
const Iterator &other)
const
1932 return !(*
this == other);
1936 template <
typename number,
bool Constness>
1938 Iterator<number, Constness>::operator<(
const Iterator &other)
const
1940 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1943 return (accessor < other.accessor);
1947 template <
typename number,
bool Constness>
1949 Iterator<number, Constness>::operator>(
const Iterator &other)
const
1951 return (other < *
this);
1955 template <
typename number,
bool Constness>
1957 Iterator<number, Constness>::operator-(
const Iterator &other)
const
1959 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1966 Iterator
copy = *
this;
1967 while (copy != other)
1975 Iterator
copy = other;
1976 while (copy != *
this)
1987 template <
typename number,
bool Constness>
1988 inline Iterator<number, Constness>
1989 Iterator<number, Constness>::operator+(
const unsigned int n)
const
1992 for (
unsigned int i = 0; i < n; ++i)
2002template <
typename number>
2006 return const_iterator(
this, 0);
2010template <
typename number>
2014 return const_iterator(
this);
2018template <
typename number>
2022 return iterator(
this, 0);
2026template <
typename number>
2030 return iterator(
this);
2034template <
typename number>
2039 return const_iterator(
this, r);
2044template <
typename number>
2049 return const_iterator(
this, r + 1);
2054template <
typename number>
2059 return iterator(
this, r);
2064template <
typename number>
2069 return iterator(
this, r + 1);
const Reference & operator=(const number n) const
const Reference & operator/=(const number n) const
const Reference & operator-=(const number n) const
const Accessor * accessor
const Reference & operator*=(const number n) const
const Reference & operator+=(const number n) const
Reference(const Accessor *accessor, const bool dummy)
MatrixType & get_matrix() const
Accessor(MatrixType *matrix)
Accessor(MatrixType *matrix, const unsigned int row)
Accessor(MatrixType *matrix, const unsigned int row)
const MatrixType & get_matrix() const
Accessor(const ChunkSparseMatrixIterators::Accessor< number, false > &a)
Accessor(MatrixType *matrix)
const ChunkSparseMatrix< number > & get_matrix() const
bool operator>(const Iterator &) const
Iterator(const ChunkSparseMatrixIterators::Iterator< number, false > &i)
const Accessor< number, Constness > & value_type
bool operator<(const Iterator &) const
Iterator(MatrixType *matrix, const unsigned int row)
const Accessor< number, Constness > * operator->() const
bool operator==(const Iterator &) const
Accessor< number, Constness > accessor
const Accessor< number, Constness > & operator*() const
int operator-(const Iterator &p) const
typename Accessor< number, Constness >::MatrixType MatrixType
bool operator!=(const Iterator &) const
Iterator operator+(const unsigned int n) const
Iterator(MatrixType *matrix)
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
somenumber matrix_scalar_product(const Vector< somenumber > &u, const Vector< somenumber > &v) const
real_type linfty_norm() const
size_type compute_location(const size_type i, const size_type j) const
ChunkSparseMatrix(const ChunkSparsityPattern &sparsity)
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
ChunkSparseMatrix< number > & operator=(const IdentityMatrix &id)
void Tvmult(OutVector &dst, const InVector &src) const
const_iterator begin() const
void add(const size_type i, const size_type j, const number value)
void set(const size_type i, const size_type j, const number value)
void precondition_SSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
const_iterator begin(const unsigned int r) const
void print(std::ostream &out) const
ChunkSparseMatrix & operator/=(const number factor)
ChunkSparseMatrix & operator=(const double d)
virtual void reinit(const ChunkSparsityPattern &sparsity)
void extract_row_copy(const size_type row, const size_type array_length, size_type &row_length, size_type *column_indices, number *values) const
iterator end(const unsigned int r)
void print_pattern(std::ostream &out, const double threshold=0.) const
void TPSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
void copy_from(const ForwardIterator begin, const ForwardIterator end)
real_type l1_norm() const
const_iterator end() const
void vmult_add(OutVector &dst, const InVector &src) const
void precondition_SOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
somenumber matrix_norm_square(const Vector< somenumber > &v) const
virtual ~ChunkSparseMatrix() override
std::unique_ptr< number[]> val
void block_write(std::ostream &out) const
void SSOR(Vector< somenumber > &v, const number omega=1.) const
void precondition_TSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
size_type n_actually_nonzero_elements() const
ChunkSparseMatrix< number > & operator=(const ChunkSparseMatrix< number > &)
ChunkSparseMatrix(const ChunkSparsityPattern &sparsity, const IdentityMatrix &id)
void vmult(OutVector &dst, const InVector &src) const
void TSOR(Vector< somenumber > &v, const number om=1.) const
void add(const number factor, const ChunkSparseMatrix< somenumber > &matrix)
size_type n_nonzero_elements() const
void copy_from(const FullMatrix< somenumber > &matrix)
ChunkSparseMatrix & operator*=(const number factor)
void block_read(std::istream &in)
number operator()(const size_type i, const size_type j) const
void Tvmult_add(OutVector &dst, const InVector &src) const
number el(const size_type i, const size_type j) const
SmartPointer< const ChunkSparsityPattern, ChunkSparseMatrix< number > > cols
ChunkSparseMatrix< number > & copy_from(const ChunkSparseMatrix< somenumber > &source)
void TSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
typename numbers::NumberTraits< number >::real_type real_type
ChunkSparseMatrix(const ChunkSparseMatrix &)
number diag_element(const size_type i) const
void SOR(Vector< somenumber > &v, const number om=1.) const
iterator begin(const unsigned int r)
void SSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
std::size_t memory_consumption() const
void add(const size_type row, const size_type n_cols, const size_type *col_indices, const number2 *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
real_type frobenius_norm() const
somenumber residual(Vector< somenumber > &dst, const Vector< somenumber > &x, const Vector< somenumber > &b) const
void PSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
void SOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
const_iterator end(const unsigned int r) const
const ChunkSparsityPattern & get_sparsity_pattern() const
Accessor(const ChunkSparsityPattern *matrix, const size_type row)
types::global_dof_index size_type
static const size_type invalid_entry
static constexpr size_type invalid_entry
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
__global__ void set(Number *val, const Number s, const size_type N)
#define DeclException0(Exception0)
static ::ExceptionBase & ExcDifferentChunkSparsityPatterns()
static ::ExceptionBase & ExcSourceEqualsDestination()
#define Assert(cond, exc)
#define AssertIsFinite(number)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcNeedsSparsityPattern()
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
static ::ExceptionBase & ExcIteratorRange(int arg1, int arg2)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
VectorType::value_type * end(VectorType &V)
void copy(const T *begin, const T *end, U *dest)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
static const bool zero_addition_can_be_elided
typename ::ChunkSparseMatrixIterators::Iterator< number, Constness >::difference_type difference_type
forward_iterator_tag iterator_category
typename ::ChunkSparseMatrixIterators::Iterator< number, Constness >::value_type value_type