16#ifndef dealii_chunk_sparse_matrix_h
17# define dealii_chunk_sparse_matrix_h
36template <
typename number>
38template <
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;
279 template <
typename,
bool>
295 template <
typename number,
bool Constness>
420template <
typename number>
679 template <
typename number2>
685 const bool elide_zero_values =
true,
686 const bool col_indices_are_sorted =
false);
731 template <
typename somenumber>
752 template <
typename ForwardIterator>
761 template <
typename somenumber>
776 template <
typename somenumber>
864 template <
class OutVector,
class InVector>
866 vmult(OutVector &dst,
const InVector &src)
const;
883 template <
class OutVector,
class InVector>
885 Tvmult(OutVector &dst,
const InVector &src)
const;
901 template <
class OutVector,
class InVector>
920 template <
class OutVector,
class InVector>
939 template <
typename somenumber>
946 template <
typename somenumber>
957 template <
typename somenumber>
1006 template <
typename somenumber>
1010 const number omega = 1.)
const;
1015 template <
typename somenumber>
1019 const number om = 1.)
const;
1024 template <
typename somenumber>
1028 const number om = 1.)
const;
1033 template <
typename somenumber>
1037 const number om = 1.)
const;
1044 template <
typename somenumber>
1052 template <
typename somenumber>
1060 template <
typename somenumber>
1074 template <
typename somenumber>
1077 const std::vector<size_type> &permutation,
1078 const std::vector<size_type> &inverse_permutation,
1079 const number om = 1.)
const;
1091 template <
typename somenumber>
1094 const std::vector<size_type> &permutation,
1095 const std::vector<size_type> &inverse_permutation,
1096 const number om = 1.)
const;
1102 template <
typename somenumber>
1106 const number om = 1.)
const;
1112 template <
typename somenumber>
1116 const number om = 1.)
const;
1122 template <
typename somenumber>
1126 const number om = 1.)
const;
1211 end(
const unsigned int r)
const;
1281 const unsigned int precision = 3,
1282 const bool scientific =
true,
1283 const unsigned int width = 0,
1284 const char * zero_string =
" ",
1285 const double denominator = 1.)
const;
1338 <<
"You are trying to access the matrix entry with index <"
1339 << arg1 <<
',' << arg2
1340 <<
">, but this entry does not exist in the sparsity pattern "
1343 "The most common cause for this problem is that you used "
1344 "a method to build the sparsity pattern that did not "
1345 "(completely) take into account all of the entries you "
1346 "will later try to write into. An example would be "
1347 "building a sparsity pattern that does not include "
1348 "the entries you will write into due to constraints "
1349 "on degrees of freedom such as hanging nodes or periodic "
1350 "boundary conditions. In such cases, building the "
1351 "sparsity pattern will succeed, but you will get errors "
1352 "such as the current one at one point or other when "
1353 "trying to write into the entries of the matrix.");
1364 <<
"The iterators denote a range of " << arg1
1365 <<
" elements, but the given number of rows was " << arg2);
1370 "You are attempting an operation on two matrices that "
1371 "are the same object, but the operation requires that the "
1372 "two objects are in fact different.");
1389 std::unique_ptr<number[]>
val;
1406 template <
typename somenumber>
1410 template <
typename,
bool>
1412 template <
typename,
bool>
1423template <
typename number>
1432template <
typename number>
1442template <
typename number>
1452template <
typename number>
1471template <
typename number>
1482 const size_type index = compute_location(i, j);
1484 ExcInvalidIndex(i, j));
1492template <
typename number>
1504 const size_type index = compute_location(i, j);
1506 ExcInvalidIndex(i, j));
1508 val[index] +=
value;
1514template <
typename number>
1515template <
typename number2>
1525 for (
size_type col = 0; col < n_cols; ++col)
1526 add(
row, col_indices[col],
static_cast<number
>(
values[col]));
1531template <
typename number>
1544 number * val_ptr = val.get();
1545 const number *
const end_ptr =
1548 while (val_ptr != end_ptr)
1549 *val_ptr++ *= factor;
1556template <
typename number>
1564 const number factor_inv = 1. / factor;
1572 number * val_ptr = val.get();
1573 const number *
const end_ptr =
1577 while (val_ptr != end_ptr)
1578 *val_ptr++ *= factor_inv;
1585template <
typename number>
1592 ExcInvalidIndex(i, j));
1593 return val[compute_location(i, j)];
1598template <
typename number>
1603 const size_type index = compute_location(i, j);
1613template <
typename number>
1631template <
typename number>
1632template <
typename ForwardIterator>
1635 const ForwardIterator
end)
1638 ExcIteratorRange(std::distance(
begin,
end), m()));
1642 using inner_iterator =
1643 typename std::iterator_traits<ForwardIterator>::value_type::const_iterator;
1645 for (ForwardIterator i =
begin; i !=
end; ++i, ++
row)
1647 const inner_iterator end_of_row = i->end();
1648 for (inner_iterator j = i->begin(); j != end_of_row; ++j)
1650 set(
row, j->first, j->second);
1661 template <
typename number>
1663 const unsigned int row)
1671 template <
typename number>
1672 inline Accessor<number, true>::Accessor(
const MatrixType *
matrix)
1679 template <
typename number>
1680 inline Accessor<number, true>::Accessor(
1683 ,
matrix(&a.get_matrix())
1688 template <
typename number>
1690 Accessor<number, true>::value()
const
1693 matrix->get_sparsity_pattern().get_chunk_size();
1700 template <
typename number>
1701 inline const typename Accessor<number, true>::MatrixType &
1702 Accessor<number, true>::get_matrix()
const
1709 template <
typename number>
1710 inline Accessor<number, false>::Reference::Reference(
const Accessor *accessor,
1712 : accessor(accessor)
1716 template <
typename number>
1717 inline Accessor<number, false>::Reference::operator number()
const
1720 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1721 return accessor->matrix
1723 accessor->chunk_row *
chunk_size + accessor->chunk_col];
1728 template <
typename number>
1729 inline const typename Accessor<number, false>::Reference &
1730 Accessor<number, false>::Reference::operator=(
const number n)
const
1733 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1736 accessor->chunk_row *
chunk_size + accessor->chunk_col] = n;
1742 template <
typename number>
1743 inline const typename Accessor<number, false>::Reference &
1744 Accessor<number, false>::Reference::operator+=(
const number n)
const
1747 accessor->matrix->get_sparsity_pattern().get_chunk_size();
1750 accessor->chunk_row *
chunk_size + accessor->chunk_col] += n;
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 Accessor<number, false>::Accessor(MatrixType *
matrix,
1800 const unsigned int row)
1808 template <
typename number>
1809 inline Accessor<number, false>::Accessor(MatrixType *
matrix)
1816 template <
typename number>
1817 inline typename Accessor<number, false>::Reference
1818 Accessor<number, false>::value()
const
1820 return Reference(
this,
true);
1825 template <
typename number>
1826 inline typename Accessor<number, false>::MatrixType &
1827 Accessor<number, false>::get_matrix()
const
1834 template <
typename number,
bool Constness>
1835 inline Iterator<number, Constness>::Iterator(MatrixType *
matrix,
1836 const unsigned int row)
1842 template <
typename number,
bool Constness>
1843 inline Iterator<number, Constness>::Iterator(MatrixType *
matrix)
1849 template <
typename number,
bool Constness>
1850 inline Iterator<number, Constness>::Iterator(
1857 template <
typename number,
bool Constness>
1858 inline Iterator<number, Constness> &
1866 template <
typename number,
bool Constness>
1867 inline Iterator<number, Constness>
1870 const Iterator iter = *
this;
1876 template <
typename number,
bool Constness>
1884 template <
typename number,
bool Constness>
1885 inline const Accessor<number, Constness> *Iterator<number, Constness>::
1892 template <
typename number,
bool Constness>
1896 return (accessor == other.accessor);
1900 template <
typename number,
bool Constness>
1904 return !(*
this == other);
1908 template <
typename number,
bool Constness>
1912 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1915 return (accessor < other.accessor);
1919 template <
typename number,
bool Constness>
1923 return (other < *
this);
1927 template <
typename number,
bool Constness>
1931 Assert(&accessor.get_matrix() == &other.accessor.get_matrix(),
1938 Iterator
copy = *
this;
1939 while (
copy != other)
1947 Iterator
copy = other;
1948 while (
copy != *
this)
1959 template <
typename number,
bool Constness>
1960 inline Iterator<number, Constness>
1964 for (
unsigned int i = 0; i < n; ++i)
1974template <
typename number>
1978 return const_iterator(
this, 0);
1982template <
typename number>
1986 return const_iterator(
this);
1990template <
typename number>
1994 return iterator(
this, 0);
1998template <
typename number>
2002 return iterator(
this);
2006template <
typename number>
2011 return const_iterator(
this, r);
2016template <
typename number>
2021 return const_iterator(
this, r + 1);
2026template <
typename number>
2031 return iterator(
this, r);
2036template <
typename number>
2041 return iterator(
this, r + 1);
bool operator!=(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
bool operator==(const AlignedVector< T > &lhs, const AlignedVector< T > &rhs)
std::enable_if< std::is_floating_point< T >::value &&std::is_floating_point< U >::value, typenameProductType< std::complex< T >, std::complex< U > >::type >::type operator*(const std::complex< T > &left, const std::complex< U > &right)
#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)
Expression operator>(const Expression &lhs, const Expression &rhs)
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
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)
number el(const size_type i, const size_type j) const
Accessor(const ChunkSparseMatrixIterators::Accessor< number, false > &a)
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
static const size_type invalid_entry
@ matrix
Contents is actually a matrix.
types::global_dof_index size_type
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(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
BarycentricPolynomial< dim, Number1 > operator-(const Number2 &a, const BarycentricPolynomial< dim, Number1 > &bp)
BarycentricPolynomial< dim, Number1 > operator+(const Number2 &a, const BarycentricPolynomial< dim, Number1 > &bp)
SynchronousIterators< Iterators > operator++(SynchronousIterators< Iterators > &a)
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)