16#ifndef dealii_chunk_sparse_matrix_h
17# define dealii_chunk_sparse_matrix_h
37template <
typename number>
39template <
typename number>
54 template <
typename number,
bool Constness>
67 template <
typename number,
bool Constness>
99 template <
typename number>
149 template <
typename,
bool>
160 template <
typename number>
198 operator number()
const;
280 template <
typename,
bool>
296 template <
typename number,
bool Constness>
415 template <
typename number,
bool Constness>
416 struct iterator_traits<
421 Iterator<number, Constness>::value_type;
423 Iterator<number, Constness>::difference_type;
447template <
typename number>
706 template <
typename number2>
711 const number2 * values,
712 const bool elide_zero_values =
true,
713 const bool col_indices_are_sorted =
false);
758 template <
typename somenumber>
779 template <
typename ForwardIterator>
788 template <
typename somenumber>
803 template <
typename somenumber>
870 number * values)
const;
891 template <
class OutVector,
class InVector>
893 vmult(OutVector &dst,
const InVector &src)
const;
910 template <
class OutVector,
class InVector>
912 Tvmult(OutVector &dst,
const InVector &src)
const;
928 template <
class OutVector,
class InVector>
947 template <
class OutVector,
class InVector>
966 template <
typename somenumber>
973 template <
typename somenumber>
984 template <
typename somenumber>
1033 template <
typename somenumber>
1037 const number omega = 1.)
const;
1042 template <
typename somenumber>
1046 const number om = 1.)
const;
1051 template <
typename somenumber>
1055 const number om = 1.)
const;
1060 template <
typename somenumber>
1064 const number om = 1.)
const;
1071 template <
typename somenumber>
1079 template <
typename somenumber>
1087 template <
typename somenumber>
1101 template <
typename somenumber>
1104 const std::vector<size_type> &permutation,
1105 const std::vector<size_type> &inverse_permutation,
1106 const number om = 1.)
const;
1118 template <
typename somenumber>
1121 const std::vector<size_type> &permutation,
1122 const std::vector<size_type> &inverse_permutation,
1123 const number om = 1.)
const;
1129 template <
typename somenumber>
1133 const number om = 1.)
const;
1139 template <
typename somenumber>
1143 const number om = 1.)
const;
1149 template <
typename somenumber>
1153 const number om = 1.)
const;
1238 end(
const unsigned int r)
const;
1308 const unsigned int precision = 3,
1309 const bool scientific =
true,
1310 const unsigned int width = 0,
1311 const char * zero_string =
" ",
1312 const double denominator = 1.)
const;
1365 <<
"You are trying to access the matrix entry with index <"
1366 << arg1 <<
',' << arg2
1367 <<
">, but this entry does not exist in the sparsity pattern "
1370 "The most common cause for this problem is that you used "
1371 "a method to build the sparsity pattern that did not "
1372 "(completely) take into account all of the entries you "
1373 "will later try to write into. An example would be "
1374 "building a sparsity pattern that does not include "
1375 "the entries you will write into due to constraints "
1376 "on degrees of freedom such as hanging nodes or periodic "
1377 "boundary conditions. In such cases, building the "
1378 "sparsity pattern will succeed, but you will get errors "
1379 "such as the current one at one point or other when "
1380 "trying to write into the entries of the matrix.");
1391 <<
"The iterators denote a range of " << arg1
1392 <<
" elements, but the given number of rows was " << arg2);
1397 "You are attempting an operation on two matrices that "
1398 "are the same object, but the operation requires that the "
1399 "two objects are in fact different.");
1416 std::unique_ptr<number[]>
val;
1433 template <
typename somenumber>
1437 template <
typename,
bool>
1439 template <
typename,
bool>
1450template <
typename number>
1459template <
typename number>
1469template <
typename number>
1479template <
typename number>
1486 cols->sparsity_pattern(i / chunk_size, j / chunk_size);
1492 return (chunk_index * chunk_size * chunk_size +
1493 (i % chunk_size) * chunk_size + (j % chunk_size));
1498template <
typename number>
1511 ExcInvalidIndex(i, j));
1519template <
typename number>
1533 ExcInvalidIndex(i, j));
1541template <
typename number>
1542template <
typename number2>
1547 const number2 * values,
1552 for (
size_type col = 0; col < n_cols; ++col)
1553 add(
row, col_indices[col],
static_cast<number
>(values[col]));
1558template <
typename number>
1571 number * val_ptr = val.get();
1572 const number *
const end_ptr =
1575 while (val_ptr != end_ptr)
1576 *val_ptr++ *= factor;
1583template <
typename number>
1591 const number factor_inv = 1. / factor;
1599 number * val_ptr = val.get();
1600 const number *
const end_ptr =
1604 while (val_ptr != end_ptr)
1605 *val_ptr++ *= factor_inv;
1612template <
typename number>
1619 ExcInvalidIndex(i, j));
1620 return val[compute_location(i, j)];
1625template <
typename number>
1640template <
typename number>
1653 (i %
chunk_size) * chunk_size + (i % chunk_size)];
1658template <
typename number>
1659template <
typename ForwardIterator>
1662 const ForwardIterator end)
1665 ExcIteratorRange(std::distance(begin, end), m()));
1669 using inner_iterator =
1670 typename std::iterator_traits<ForwardIterator>::value_type::const_iterator;
1672 for (ForwardIterator i = begin; i !=
end; ++i, ++
row)
1674 const inner_iterator end_of_row = i->end();
1675 for (inner_iterator j = i->begin(); j != end_of_row; ++j)
1677 set(
row, j->first, j->second);
1688 template <
typename number>
1690 const unsigned int row)
1698 template <
typename number>
1699 inline Accessor<number, true>::Accessor(
const MatrixType *matrix)
1706 template <
typename number>
1707 inline Accessor<number, true>::Accessor(
1710 ,
matrix(&a.get_matrix())
1715 template <
typename number>
1717 Accessor<number, true>::value()
const
1720 matrix->get_sparsity_pattern().get_chunk_size();
1727 template <
typename number>
1728 inline const typename Accessor<number, true>::MatrixType &
1729 Accessor<number, true>::get_matrix()
const
1736 template <
typename number>
1737 inline Accessor<number, false>::Reference::Reference(
const Accessor *accessor,
1739 : accessor(accessor)
1743 template <
typename number>
1744 inline Accessor<number, false>::Reference::operator number()
const
1747 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1748 return accessor->matrix
1750 accessor->chunk_row *
chunk_size + accessor->chunk_col];
1755 template <
typename number>
1756 inline const typename Accessor<number, false>::Reference &
1757 Accessor<number, false>::Reference::operator=(
const number n)
const
1760 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1763 accessor->chunk_row *
chunk_size + accessor->chunk_col] = n;
1769 template <
typename number>
1770 inline const typename Accessor<number, false>::Reference &
1771 Accessor<number, false>::Reference::operator+=(
const number n)
const
1774 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1777 accessor->chunk_row *
chunk_size + accessor->chunk_col] += n;
1783 template <
typename number>
1784 inline const typename Accessor<number, false>::Reference &
1785 Accessor<number, false>::Reference::operator-=(
const number n)
const
1788 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1791 accessor->chunk_row *
chunk_size + accessor->chunk_col] -= n;
1797 template <
typename number>
1798 inline const typename Accessor<number, false>::Reference &
1799 Accessor<number, false>::Reference::operator*=(
const number n)
const
1802 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1805 accessor->chunk_row *
chunk_size + accessor->chunk_col] *= n;
1811 template <
typename number>
1812 inline const typename Accessor<number, false>::Reference &
1813 Accessor<number, false>::Reference::operator/=(
const number n)
const
1816 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1819 accessor->chunk_row *
chunk_size + accessor->chunk_col] /= n;
1825 template <
typename number>
1826 inline Accessor<number, false>::Accessor(MatrixType * matrix,
1827 const unsigned int row)
1835 template <
typename number>
1836 inline Accessor<number, false>::Accessor(MatrixType *matrix)
1843 template <
typename number>
1844 inline typename Accessor<number, false>::Reference
1845 Accessor<number, false>::value()
const
1847 return Reference(
this,
true);
1852 template <
typename number>
1853 inline typename Accessor<number, false>::MatrixType &
1854 Accessor<number, false>::get_matrix()
const
1861 template <
typename number,
bool Constness>
1862 inline Iterator<number, Constness>::Iterator(MatrixType * matrix,
1863 const unsigned int row)
1869 template <
typename number,
bool Constness>
1870 inline Iterator<number, Constness>::Iterator(MatrixType *matrix)
1876 template <
typename number,
bool Constness>
1877 inline Iterator<number, Constness>::Iterator(
1884 template <
typename number,
bool Constness>
1885 inline Iterator<number, Constness> &
1886 Iterator<number, Constness>::operator++()
1893 template <
typename number,
bool Constness>
1894 inline Iterator<number, Constness>
1895 Iterator<number, Constness>::operator++(
int)
1897 const Iterator iter = *
this;
1903 template <
typename number,
bool Constness>
1904 inline const Accessor<number, Constness> &
1905 Iterator<number, Constness>::operator*()
const
1911 template <
typename number,
bool Constness>
1912 inline const Accessor<number, Constness> *
1913 Iterator<number, Constness>::operator->()
const
1919 template <
typename number,
bool Constness>
1921 Iterator<number, Constness>::operator==(
const Iterator &other)
const
1923 return (accessor == other.accessor);
1927 template <
typename number,
bool Constness>
1929 Iterator<number, Constness>::operator!=(
const Iterator &other)
const
1931 return !(*
this == other);
1935 template <
typename number,
bool Constness>
1937 Iterator<number, Constness>::operator<(
const Iterator &other)
const
1939 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1942 return (accessor < other.accessor);
1946 template <
typename number,
bool Constness>
1948 Iterator<number, Constness>::operator>(
const Iterator &other)
const
1950 return (other < *
this);
1954 template <
typename number,
bool Constness>
1956 Iterator<number, Constness>::operator-(
const Iterator &other)
const
1958 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1965 Iterator
copy = *
this;
1966 while (copy != other)
1974 Iterator
copy = other;
1975 while (copy != *
this)
1986 template <
typename number,
bool Constness>
1987 inline Iterator<number, Constness>
1988 Iterator<number, Constness>::operator+(
const unsigned int n)
const
1991 for (
unsigned int i = 0; i < n; ++i)
2001template <
typename number>
2005 return const_iterator(
this, 0);
2009template <
typename number>
2013 return const_iterator(
this);
2017template <
typename number>
2021 return iterator(
this, 0);
2025template <
typename number>
2029 return iterator(
this);
2033template <
typename number>
2038 return const_iterator(
this, r);
2043template <
typename number>
2048 return const_iterator(
this, r + 1);
2053template <
typename number>
2058 return iterator(
this, r);
2063template <
typename number>
2068 return iterator(
this, r + 1);
static const 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)
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 DeclException0(Exception0)
size_type compute_location(const size_type i, const size_type j) const
static ::ExceptionBase & ExcDifferentChunkSparsityPatterns()
void print(std::ostream &out) const
static ::ExceptionBase & ExcSourceEqualsDestination()
void print_pattern(std::ostream &out, const double threshold=0.) const
#define Assert(cond, exc)
std::unique_ptr< number[]> val
void block_write(std::ostream &out) const
#define AssertIsFinite(number)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcNeedsSparsityPattern()
void block_read(std::istream &in)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
SmartPointer< const ChunkSparsityPattern, ChunkSparseMatrix< number > > cols
static ::ExceptionBase & ExcNotInitialized()
static ::ExceptionBase & ExcNotQuadratic()
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
static ::ExceptionBase & ExcIteratorRange(int arg1, int arg2)
static const bool zero_addition_can_be_elided
#define AssertThrow(cond, exc)
bool operator>(const Iterator &) const
somenumber matrix_scalar_product(const Vector< somenumber > &u, const Vector< somenumber > &v) const
Iterator(const ChunkSparseMatrixIterators::Iterator< number, false > &i)
Accessor(MatrixType *matrix, const unsigned int row)
const Accessor< number, Constness > & value_type
real_type linfty_norm() const
ChunkSparseMatrix(const ChunkSparsityPattern &sparsity)
bool operator<(const Iterator &) const
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
Iterator(MatrixType *matrix, const unsigned int row)
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
const Accessor< number, Constness > * operator->() const
MatrixType & get_matrix() const
ChunkSparseMatrix & operator/=(const number factor)
ChunkSparseMatrix & operator=(const double d)
virtual void reinit(const ChunkSparsityPattern &sparsity)
const MatrixType & get_matrix() 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
iterator end(const unsigned int r)
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
typename ::ChunkSparseMatrixIterators::Iterator< number, Constness >::difference_type difference_type
bool operator==(const Iterator &) 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
Accessor< number, Constness > accessor
somenumber matrix_norm_square(const Vector< somenumber > &v) const
Accessor(MatrixType *matrix)
virtual ~ChunkSparseMatrix() override
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
const Reference & operator=(const number n) const
const Accessor< number, Constness > & operator*() const
const Reference & operator/=(const number n) const
const Reference & operator-=(const number n) 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
int operator-(const Iterator &p) const
const Accessor * accessor
const ChunkSparseMatrix< number > & get_matrix() const
void add(const number factor, const ChunkSparseMatrix< somenumber > &matrix)
size_type n_nonzero_elements() const
void copy_from(const FullMatrix< somenumber > &matrix)
const Reference & operator*=(const number n) const
ChunkSparseMatrix & operator*=(const number factor)
number operator()(const size_type i, const size_type j) const
void Tvmult_add(OutVector &dst, const InVector &src) const
const Reference & operator+=(const number n) const
Reference(const Accessor *accessor, const bool dummy)
forward_iterator_tag iterator_category
number el(const size_type i, const size_type j) const
Accessor(const ChunkSparseMatrixIterators::Accessor< number, false > &a)
typename ::ChunkSparseMatrixIterators::Iterator< number, Constness >::value_type value_type
Accessor(MatrixType *matrix)
types::global_dof_index size_type
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
typename Accessor< number, Constness >::MatrixType MatrixType
Accessor(MatrixType *matrix, const unsigned int row)
std::size_t memory_consumption() const
bool operator!=(const Iterator &) 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)
Iterator operator+(const unsigned int n) const
Iterator(MatrixType *matrix)
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)
static const size_type invalid_entry
types::global_dof_index size_type
@ 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