15#ifndef dealii_trilinos_tpetra_sparse_matrix_h
16#define dealii_trilinos_tpetra_sparse_matrix_h
20#ifdef DEAL_II_TRILINOS_WITH_TPETRA
33# include <Tpetra_Core.hpp>
34# include <Tpetra_CrsMatrix.hpp>
36# include <type_traits>
43template <
typename Number>
47template <
typename MatrixType>
52 namespace TpetraWrappers
54 template <
typename MemorySpace>
59 template <
typename Number,
typename MemorySpace,
bool Constness>
69 namespace TpetraWrappers
109 template <
typename Number,
typename MemorySpace = ::MemorySpace::Host>
123 <<
"You tried to access row " << arg1
124 <<
" of a non-contiguous locally owned row set."
125 <<
" The row " << arg1
126 <<
" is not stored locally and can't be accessed.");
185 const unsigned int n_max_entries_per_row);
196 const std::vector<unsigned int> &n_entries_per_row);
242 template <
typename SparsityPatternType>
244 reinit(
const SparsityPatternType &sparsity_pattern);
276 const MPI_Comm communicator = MPI_COMM_WORLD,
277 const unsigned int n_max_entries_per_row = 0);
289 const std::vector<unsigned int> &n_entries_per_row);
306 const IndexSet &col_parallel_partitioning,
307 const MPI_Comm communicator = MPI_COMM_WORLD,
308 const size_type n_max_entries_per_row = 0);
319 const IndexSet &col_parallel_partitioning,
321 const std::vector<unsigned int> &n_entries_per_row);
342 template <
typename SparsityPatternType>
344 !std::is_same_v<SparsityPatternType, ::SparseMatrix<double>>>
346 const SparsityPatternType &sparsity_pattern,
347 const MPI_Comm communicator = MPI_COMM_WORLD,
348 const bool exchange_data =
false);
362 template <
typename SparsityPatternType>
364 !std::is_same_v<SparsityPatternType, ::SparseMatrix<double>>>
366 const IndexSet &col_parallel_partitioning,
367 const SparsityPatternType &sparsity_pattern,
368 const MPI_Comm communicator = MPI_COMM_WORLD,
369 const bool exchange_data =
false);
389 const IndexSet &col_parallel_partitioning,
390 const ::SparseMatrix<Number> &dealii_sparse_matrix,
391 const MPI_Comm communicator = MPI_COMM_WORLD,
392 const double drop_tolerance = 1e-13,
393 const bool copy_values =
true,
394 const ::SparsityPattern *use_this_sparsity =
nullptr);
434 std::pair<size_type, size_type>
536 const Number *values,
537 const bool elide_zero_values =
true,
538 const bool col_indices_are_sorted =
false);
607 set(
const std::vector<size_type> &indices,
609 const bool elide_zero_values =
false);
617 set(
const std::vector<size_type> &row_indices,
618 const std::vector<size_type> &col_indices,
620 const bool elide_zero_values =
false);
651 const std::vector<size_type> &col_indices,
652 const std::vector<Number> &values,
653 const bool elide_zero_values =
false);
686 const Number *values,
687 const bool elide_zero_values =
false);
740 const Number new_diag_value = 0);
812 template <
typename InputVectorType>
814 vmult(InputVectorType &dst,
const InputVectorType &src)
const;
823 template <
typename InputVectorType>
825 Tvmult(InputVectorType &dst,
const InputVectorType &src)
const;
833 template <
typename InputVectorType>
835 vmult_add(InputVectorType &dst,
const InputVectorType &src)
const;
844 template <
typename InputVectorType>
846 Tvmult_add(InputVectorType &dst,
const InputVectorType &src)
const;
945 const bool print_detailed_trilinos_information =
false)
const;
1017 Teuchos::RCP<const TpetraTypes::MatrixType<Number, MemorySpace>>
1029 Teuchos::RCP<TpetraTypes::MatrixType<Number, MemorySpace>>
1171 "You are attempting an operation on two vectors that "
1172 "are the same object, but the operation requires that the "
1173 "two objects are in fact different.");
1179 "The column partitioning of a matrix does not match "
1180 "the partitioning of a vector you are trying to "
1181 "multiply it with. Are you multiplying the "
1182 "matrix with a vector that has ghost elements?");
1188 "The row partitioning of a matrix does not match "
1189 "the partitioning of a vector you are trying to "
1190 "put the result of a matrix-vector product in. "
1191 "Are you trying to put the product of the "
1192 "matrix with a vector into a vector that has "
1201 <<
"The entry with index <" << arg1 <<
',' << arg2
1202 <<
"> does not exist.");
1212 <<
"You tried to access element (" << arg1 <<
'/' << arg2
1214 <<
" of a distributed matrix, but only rows in range ["
1215 << arg3 <<
',' << arg4
1216 <<
"] are stored locally and can be accessed.");
1243 Teuchos::RCP<TpetraTypes::MatrixType<Number, MemorySpace>>
matrix;
1297 template <
typename Number,
typename MemorySpace>
1369 std::shared_ptr<std::vector<::types::signed_global_dof_index>>
1392 template <
typename Number,
typename MemorySpace,
bool Constness>
1411 template <
typename Number,
typename MemorySpace>
1439 template <
bool Other>
1452 template <
typename Number,
typename MemorySpace>
1467 operator Number()
const;
1473 operator=(
const Number n)
const;
1479 operator+=(
const Number n)
const;
1485 operator-=(
const Number n)
const;
1491 operator*=(
const Number n)
const;
1497 operator/=(
const Number n)
const;
1551 template <
typename Number,
typename MemorySpace,
bool Constness>
1590 template <
bool Other>
1621 template <
bool OtherConstness>
1628 template <
bool OtherConstness>
1637 template <
bool OtherConstness>
1644 template <
bool OtherConstness>
1654 <<
"Attempt to access element " << arg2 <<
" of row "
1655 << arg1 <<
" which doesn't have that many elements.");
1677 template <
typename Number,
typename MemorySpace,
bool Constness>
1680 Iterator<Number, MemorySpace, Constness>>
1684 typename ::LinearAlgebra::TpetraWrappers::SparseMatrixIterators::
1685 Iterator<Number, MemorySpace, Constness>::value_type;
1687 typename ::LinearAlgebra::TpetraWrappers::SparseMatrixIterators::
1688 Iterator<Number, MemorySpace, Constness>::difference_type;
1700 namespace TpetraWrappers
1702 template <
typename Number,
typename MemorySpace>
1708 set(i, 1, &j, &value,
false);
1713 template <
typename Number,
typename MemorySpace>
1719 add(i, 1, &j, &value,
false);
1724 template <
typename Number,
typename MemorySpace>
1740 template <
typename Number,
typename MemorySpace>
1744 return matrix->getRowMap()->getGlobalNumElements();
1749 template <
typename Number,
typename MemorySpace>
1757 return column_space_map->getGlobalNumElements();
1762 template <
typename Number,
typename MemorySpace>
1771 template <
typename Number,
typename MemorySpace>
1780 template <
typename Number,
typename MemorySpace>
1789 template <
typename Number,
typename MemorySpace>
1794 if (in_local_range(r) && (row_length(r) > 0))
1801 template <
typename Number,
typename MemorySpace>
1811 if (in_local_range(i) && (row_length(i) > 0))
1821 template <
typename Number,
typename MemorySpace>
1830 template <
typename Number,
typename MemorySpace>
1839 template <
typename Number,
typename MemorySpace>
1844 if (in_local_range(r) && (row_length(r) > 0))
1852 template <
typename Number,
typename MemorySpace>
1862 if (in_local_range(i) && (row_length(i) > 0))
1872 template <
typename Number,
typename MemorySpace>
1877 const size_type begin = matrix->getRowMap()->getMinGlobalIndex();
1878 const size_type end = matrix->getRowMap()->getMaxGlobalIndex() + 1;
1880 return ((index >= begin) && (index < end));
1885 template <
typename Number,
typename MemorySpace>
1889 auto n_entries = matrix->getNumEntriesInGlobalRow(row);
1891 Teuchos::OrdinalTraits<
decltype(n_entries)>::invalid(),
1892 ExcAccessToNonlocalRow(row));
1899 template <
typename Number,
typename MemorySpace>
1908 template <
typename Number,
typename MemorySpace>
1917 template <
typename Number,
typename MemorySpace>
1926 template <
typename Number,
typename MemorySpace>
1935 template <
typename Number,
typename MemorySpace>
1936 inline Teuchos::RCP<const TpetraTypes::MatrixType<Number, MemorySpace>>
1939 return matrix.getConst();
1944 template <
typename Number,
typename MemorySpace>
1945 inline Teuchos::RCP<TpetraTypes::MatrixType<Number, MemorySpace>>
1953 template <
typename Number,
typename MemorySpace>
1957 return IndexSet(matrix->getDomainMap());
1962 template <
typename Number,
typename MemorySpace>
1966 return IndexSet(matrix->getRangeMap());
1972 template <
typename Number,
typename MemorySpace>
1985 template <
typename Number,
typename MemorySpace>
1994 template <
typename Number,
typename MemorySpace>
1999 return (*colnum_cache)[a_index];
2003 template <
typename Number,
typename MemorySpace>
2013 template <
typename Number,
typename MemorySpace>
2023 if ((this->a_row == matrix->m()) ||
2024 (matrix->in_local_range(this->a_row) ==
false))
2026 colnum_cache.reset();
2027 value_cache.reset();
2035 matrix->row_length(this->a_row);
2036 if (value_cache.get() ==
nullptr)
2038 value_cache = std::make_shared<std::vector<Number>>(colnums);
2039 colnum_cache = std::make_shared<
2040 std::vector<::types::signed_global_dof_index>>(colnums);
2044 value_cache->resize(colnums);
2045 colnum_cache->resize(colnums);
2048# if DEAL_II_TRILINOS_VERSION_GTE(13, 2, 0)
2050 nonconst_global_inds_host_view_type col_indices(colnum_cache->data(),
2053 nonconst_values_host_view_type values(value_cache->data(), colnums);
2055 Teuchos::ArrayView<::types::signed_global_dof_index> col_indices(
2057 Teuchos::ArrayView<Number> values(*value_cache);
2059 matrix->trilinos_matrix().getGlobalRowCopy(this->a_row,
2076 template <
typename Number,
typename MemorySpace>
2089 template <
typename Number,
typename MemorySpace>
2090 template <
bool Other>
2098 template <
typename Number,
typename MemorySpace>
2111 template <
typename Number,
typename MemorySpace>
2119 template <
typename Number,
typename MemorySpace>
2123 return (*accessor.value_cache)[accessor.a_index];
2128 template <
typename Number,
typename MemorySpace>
2131 const Number n)
const
2133 (*accessor.value_cache)[accessor.a_index] = n;
2134 accessor.matrix->set(accessor.row(),
2136 static_cast<Number
>(*
this));
2142 template <
typename Number,
typename MemorySpace>
2145 const Number n)
const
2147 (*accessor.value_cache)[accessor.a_index] += n;
2148 accessor.matrix->set(accessor.row(),
2150 static_cast<Number
>(*
this));
2156 template <
typename Number,
typename MemorySpace>
2159 const Number n)
const
2161 (*accessor.value_cache)[accessor.a_index] -= n;
2162 accessor.matrix->set(accessor.row(),
2164 static_cast<Number
>(*
this));
2169 template <
typename Number,
typename MemorySpace>
2172 const Number n)
const
2174 (*accessor.value_cache)[accessor.a_index] *= n;
2175 accessor.matrix->set(accessor.row(),
2177 static_cast<Number
>(*
this));
2182 template <
typename Number,
typename MemorySpace>
2185 const Number n)
const
2187 (*accessor.value_cache)[accessor.a_index] /= n;
2188 accessor.matrix->set(accessor.row(),
2190 static_cast<Number
>(*
this));
2195 template <
typename Number,
typename MemorySpace>
2204 template <
typename Number,
typename MemorySpace>
2216 template <
typename Number,
typename MemorySpace,
bool Constness>
2221 : accessor(matrix, row, index)
2225 template <
typename Number,
typename MemorySpace,
bool Constness>
2226 template <
bool Other>
2229 : accessor(other.accessor)
2234 template <
typename Number,
typename MemorySpace,
bool Constness>
2246 if (accessor.a_index >= accessor.colnum_cache->size())
2248 accessor.a_index = 0;
2252 (accessor.a_row < accessor.matrix->m()) &&
2253 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2254 (accessor.matrix->row_length(accessor.a_row) == 0)))
2257 accessor.visit_present_row();
2263 template <
typename Number,
typename MemorySpace,
bool Constness>
2274 template <
typename Number,
typename MemorySpace,
bool Constness>
2283 template <
typename Number,
typename MemorySpace,
bool Constness>
2292 template <
typename Number,
typename MemorySpace,
bool Constness>
2293 template <
bool OtherConstness>
2298 return (accessor.a_row == other.
accessor.a_row &&
2299 accessor.a_index == other.
accessor.a_index);
2304 template <
typename Number,
typename MemorySpace,
bool Constness>
2305 template <
bool OtherConstness>
2310 return !(*
this == other);
2315 template <
typename Number,
typename MemorySpace,
bool Constness>
2316 template <
bool OtherConstness>
2321 return (accessor.row() < other.
accessor.row() ||
2322 (accessor.row() == other.
accessor.row() &&
2323 accessor.index() < other.
accessor.index()));
2327 template <
typename Number,
typename MemorySpace,
bool Constness>
2328 template <
bool OtherConstness>
2333 return (other < *
this);
SparseMatrix< Number, MemorySpace > * matrix
std::shared_ptr< std::vector< Number > > value_cache
std::shared_ptr< std::vector<::types::signed_global_dof_index > > colnum_cache
AccessorBase(SparseMatrix< Number, MemorySpace > *matrix, const size_type row, const size_type index)
typename AccessorBase< Number, MemorySpace >::size_type size_type
typename AccessorBase< Number, MemorySpace >::size_type size_type
bool operator!=(const Iterator< Number, MemorySpace, OtherConstness > &) const
bool operator<(const Iterator< Number, MemorySpace, OtherConstness > &) const
bool operator>(const Iterator< Number, MemorySpace, OtherConstness > &) const
::types::global_dof_index size_type
typename Accessor< Number, MemorySpace, Constness >::MatrixType MatrixType
Iterator(MatrixType *matrix, const size_type row, const size_type index)
Accessor< Number, MemorySpace, Constness > accessor
const Accessor< Number, MemorySpace, Constness > * operator->() const
Iterator< Number, MemorySpace, Constness > & operator++()
const Accessor< Number, MemorySpace, Constness > & operator*() const
bool operator==(const Iterator< Number, MemorySpace, OtherConstness > &) const
virtual ~SparseMatrix() override=default
unsigned int row_length(const size_type row) const
SparseMatrix(const SparseMatrix< Number, MemorySpace > &)=delete
void clear_row(const size_type row, const Number new_diag_value=0)
SparseMatrix(const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const size_type n_max_entries_per_row=0)
void print(std::ostream &out, const bool print_detailed_trilinos_information=false) const
IndexSet locally_owned_domain_indices() const
::types::global_dof_index size_type
iterator begin(const size_type r)
SparseMatrix< Number, MemorySpace > & operator=(SparseMatrix< Number, MemorySpace > &&other) noexcept
TpetraTypes::MatrixType< Number, MemorySpace > & trilinos_matrix()
void add(const size_type i, const size_type j, const Number value)
void set(const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, const FullMatrix< Number > &full_matrix, const bool elide_zero_values=false)
SparseMatrix(SparseMatrix< Number, MemorySpace > &&other) noexcept
void clear_rows(const ArrayView< const size_type > &rows, const Number new_diag_value=0)
const_iterator end() const
void set(const std::vector< size_type > &indices, const FullMatrix< Number > &full_matrix, const bool elide_zero_values=false)
Number residual(Vector< Number, MemorySpace > &dst, const Vector< Number, MemorySpace > &x, const Vector< Number, MemorySpace > &b) const
std::pair< size_type, size_type > local_range() const
void set(const size_type row, const std::vector< size_type > &col_indices, const std::vector< Number > &values, const bool elide_zero_values=false)
Number frobenius_norm() const
void Tvmult(InputVectorType &dst, const InputVectorType &src) const
Number matrix_scalar_product(const Vector< Number, MemorySpace > &u, const Vector< Number, MemorySpace > &v) const
Teuchos::RCP< const TpetraTypes::MatrixType< Number, MemorySpace > > trilinos_rcp() const
Teuchos::RCP< TpetraTypes::MatrixType< Number, MemorySpace > > matrix
void reinit(const SparsityPatternType &sparsity_pattern)
SparseMatrix & operator=(const double d)
Number element(const size_type i, const size_type j, const bool no_error) const
SparseMatrix(const size_type m, const size_type n, const unsigned int n_max_entries_per_row)
SparseMatrix & operator/=(const Number factor)
Teuchos::RCP< TpetraTypes::MapType< MemorySpace > > column_space_map
SparseMatrix & operator*=(const Number factor)
bool in_local_range(const size_type index) const
MPI_Comm get_mpi_communicator() const
void add(const size_type row, const size_type n_cols, const size_type *col_indices, const Number *values, const bool elide_zero_values=true, const bool col_indices_are_sorted=false)
const_iterator begin(const size_type r) const
const_iterator begin() const
SparseMatrix(const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const MPI_Comm communicator, const std::vector< unsigned int > &n_entries_per_row)
SparseMatrix(const SparsityPattern< MemorySpace > &sparsity_pattern)
SparseMatrix(const IndexSet ¶llel_partitioning, const MPI_Comm communicator=MPI_COMM_WORLD, const unsigned int n_max_entries_per_row=0)
Number diag_element(const size_type i) const
const_iterator end(const size_type r) const
Number operator()(const size_type i, const size_type j) const
void vmult(InputVectorType &dst, const InputVectorType &src) const
SparseMatrix(const size_type m, const size_type n, const std::vector< unsigned int > &n_entries_per_row)
void set(const size_type i, const size_type j, const Number value)
void copy_from(const SparseMatrix< Number, MemorySpace > &source)
iterator end(const size_type r)
SparseMatrix(const IndexSet ¶llel_partitioning, const MPI_Comm communicator, const std::vector< unsigned int > &n_entries_per_row)
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > reinit(const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
SparseMatrix< Number, MemorySpace > & operator=(const SparseMatrix< Number, MemorySpace > &)=delete
IndexSet locally_owned_range_indices() const
Teuchos::RCP< TpetraTypes::MatrixType< Number, MemorySpace > > trilinos_rcp()
size_t n_nonzero_elements() const
void Tvmult_add(InputVectorType &dst, const InputVectorType &src) const
void reinit(const SparsityPattern< MemorySpace > &sparsity_pattern)
Number el(const size_type i, const size_type j) const
void vmult_add(InputVectorType &dst, const InputVectorType &src) const
const TpetraTypes::MatrixType< Number, MemorySpace > & trilinos_matrix() const
void compress(VectorOperation::values operation)
std::enable_if_t< !std::is_same_v< SparsityPatternType, ::SparseMatrix< double > > > reinit(const IndexSet ¶llel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
void add(const Number factor, const SparseMatrix< Number, MemorySpace > &matrix)
void reinit(const IndexSet &row_parallel_partitioning, const IndexSet &col_parallel_partitioning, const ::SparseMatrix< Number > &dealii_sparse_matrix, const MPI_Comm communicator=MPI_COMM_WORLD, const double drop_tolerance=1e-13, const bool copy_values=true, const ::SparsityPattern *use_this_sparsity=nullptr)
void set(const size_type row, const size_type n_cols, const size_type *col_indices, const Number *values, const bool elide_zero_values=false)
bool is_compressed() const
Number matrix_norm_square(const Vector< Number, MemorySpace > &v) const
unsigned int local_size() const
real_type l2_norm() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcDomainMapMismatch()
static ::ExceptionBase & ExcColMapMismatch()
#define DeclException0(Exception0)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
#define Assert(cond, exc)
static ::ExceptionBase & ExcMatrixNotCompressed()
static ::ExceptionBase & ExcIteratorPastEnd()
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcBeyondEndOfMatrix()
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcSourceEqualsDestination()
Tpetra::CrsMatrix< Number, LO, GO, NodeType< MemorySpace > > MatrixType
unsigned int global_dof_index
static const bool zero_addition_can_be_elided
forward_iterator_tag iterator_category
typename ::LinearAlgebra::TpetraWrappers::SparseMatrixIterators:: Iterator< Number, MemorySpace, Constness >::value_type value_type
typename ::LinearAlgebra::TpetraWrappers::SparseMatrixIterators:: Iterator< Number, MemorySpace, Constness >::difference_type difference_type