16#ifndef dealii_trilinos_sparse_matrix_h
17# define dealii_trilinos_sparse_matrix_h
22# ifdef DEAL_II_WITH_TRILINOS
37# include <Epetra_Comm.h>
38# include <Epetra_CrsGraph.h>
39# include <Epetra_Export.h>
40# include <Epetra_FECrsMatrix.h>
41# include <Epetra_Map.h>
42# include <Epetra_MpiComm.h>
43# include <Epetra_MultiVector.h>
44# include <Epetra_Operator.h>
49# include <type_traits>
56template <
typename MatrixType>
59template <
typename number>
71 template <
bool Constness>
96 <<
"You tried to access row " << arg1
97 <<
" of a distributed sparsity pattern, "
98 <<
" but only rows " << arg2 <<
" through " << arg3
99 <<
" are stored locally and can be accessed.");
197 template <
bool Constess>
238 template <
bool Other>
335 friend class Reference;
351 template <
bool Constness>
387 template <
bool Other>
418 template <
bool OtherConstness>
425 template <
bool OtherConstness>
434 template <
bool OtherConstness>
441 template <
bool OtherConstness>
451 <<
"Attempt to access element " << arg2 <<
" of row "
452 << arg1 <<
" which doesn't have that many elements.");
460 template <
bool Other>
471 template <
bool Constness>
477 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
478 Constness>::value_type;
480 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
481 Constness>::difference_type;
563 <<
"You tried to access row " << arg1
564 <<
" of a non-contiguous locally owned row set."
565 <<
" The row " << arg1
566 <<
" is not stored locally and can't be accessed.");
617 const unsigned int n_max_entries_per_row);
628 const std::vector<unsigned int> &n_entries_per_row);
672 template <
typename SparsityPatternType>
674 reinit(
const SparsityPatternType &sparsity_pattern);
722 template <
typename number>
724 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
725 const double drop_tolerance = 1e-13,
726 const bool copy_values =
true,
727 const ::SparsityPattern * use_this_sparsity =
nullptr);
735 reinit(
const Epetra_CrsMatrix &input_matrix,
const bool copy_values =
true);
755 const MPI_Comm communicator = MPI_COMM_WORLD,
756 const unsigned int n_max_entries_per_row = 0);
767 const std::vector<unsigned int> &n_entries_per_row);
784 const IndexSet &col_parallel_partitioning,
785 const MPI_Comm communicator = MPI_COMM_WORLD,
786 const size_type n_max_entries_per_row = 0);
803 const IndexSet & col_parallel_partitioning,
805 const std::vector<unsigned int> &n_entries_per_row);
827 template <
typename SparsityPatternType>
830 const SparsityPatternType &sparsity_pattern,
831 const MPI_Comm communicator = MPI_COMM_WORLD,
832 const bool exchange_data =
false);
846 template <
typename SparsityPatternType>
848 !std::is_same<SparsityPatternType, ::SparseMatrix<double>>::value>
850 const IndexSet & col_parallel_partitioning,
851 const SparsityPatternType &sparsity_pattern,
852 const MPI_Comm communicator = MPI_COMM_WORLD,
853 const bool exchange_data =
false);
871 template <
typename number>
874 const ::SparseMatrix<number> &dealii_sparse_matrix,
875 const MPI_Comm communicator = MPI_COMM_WORLD,
876 const double drop_tolerance = 1e-13,
877 const bool copy_values =
true,
878 const ::SparsityPattern * use_this_sparsity =
nullptr);
893 template <
typename number>
896 const IndexSet & col_parallel_partitioning,
897 const ::SparseMatrix<number> &dealii_sparse_matrix,
898 const MPI_Comm communicator = MPI_COMM_WORLD,
899 const double drop_tolerance = 1e-13,
900 const bool copy_values =
true,
901 const ::SparsityPattern * use_this_sparsity =
nullptr);
939 std::pair<size_type, size_type>
1100 set(
const std::vector<size_type> & indices,
1102 const bool elide_zero_values =
false);
1110 set(
const std::vector<size_type> & row_indices,
1111 const std::vector<size_type> & col_indices,
1113 const bool elide_zero_values =
false);
1144 const std::vector<size_type> & col_indices,
1145 const std::vector<TrilinosScalar> &values,
1146 const bool elide_zero_values =
false);
1175 template <
typename Number>
1180 const Number * values,
1181 const bool elide_zero_values =
false);
1214 add(
const std::vector<size_type> & indices,
1216 const bool elide_zero_values =
true);
1224 add(
const std::vector<size_type> & row_indices,
1225 const std::vector<size_type> & col_indices,
1227 const bool elide_zero_values =
true);
1244 const std::vector<size_type> & col_indices,
1245 const std::vector<TrilinosScalar> &values,
1246 const bool elide_zero_values =
true);
1266 const bool elide_zero_values =
true,
1267 const bool col_indices_are_sorted =
false);
1347 clear_rows(
const std::vector<size_type> &rows,
1443 template <
typename VectorType>
1445 std::is_same<typename VectorType::value_type, TrilinosScalar>::value>
1446 vmult(VectorType &dst,
const VectorType &src)
const;
1454 template <
typename VectorType>
1456 !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>
1457 vmult(VectorType &dst,
const VectorType &src)
const;
1473 template <
typename VectorType>
1475 std::is_same<typename VectorType::value_type, TrilinosScalar>::value>
1476 Tvmult(VectorType &dst,
const VectorType &src)
const;
1484 template <
typename VectorType>
1486 !std::is_same<typename VectorType::value_type, TrilinosScalar>::value>
1487 Tvmult(VectorType &dst,
const VectorType &src)
const;
1498 template <
typename VectorType>
1500 vmult_add(VectorType &dst,
const VectorType &src)
const;
1512 template <
typename VectorType>
1514 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
1668 const Epetra_CrsMatrix &
1675 const Epetra_CrsGraph &
1824 print(std::ostream &out,
1825 const bool write_extended_trilinos_info =
false)
const;
1838 <<
"An error with error number " << arg1
1839 <<
" occurred while calling a Trilinos function. "
1841 "For historical reasons, many Trilinos functions express "
1842 "errors by returning specific integer values to indicate "
1843 "certain errors. Unfortunately, different Trilinos functions "
1844 "often use the same integer values for different kinds of "
1845 "errors, and in most cases it is also not documented what "
1846 "each error code actually means. As a consequence, it is often "
1847 "difficult to say what a particular error (in this case, "
1848 "the error with integer code '"
1850 <<
"') represents and how one should fix a code to avoid it. "
1851 "The best one can often do is to look up the call stack to "
1852 "see which deal.II function generated the error, and which "
1853 "Trilinos function the error code had originated from; "
1854 "then look up the Trilinos source code of that function (for "
1855 "example on github) to see what code path set that error "
1856 "code. Short of going through all of that, the only other "
1857 "option is to guess the cause of the error from "
1858 "the context in which the error appeared.");
1867 <<
"The entry with index <" << arg1 <<
',' << arg2
1868 <<
"> does not exist.");
1874 "You are attempting an operation on two vectors that "
1875 "are the same object, but the operation requires that the "
1876 "two objects are in fact different.");
1891 <<
"You tried to access element (" << arg1 <<
'/' << arg2
1893 <<
" of a distributed matrix, but only rows in range ["
1894 << arg3 <<
',' << arg4
1895 <<
"] are stored locally and can be accessed.");
1903 <<
"You tried to access element (" << arg1 <<
'/' << arg2
1904 <<
')' <<
" of a sparse matrix, but it appears to not"
1905 <<
" exist in the Trilinos sparsity pattern.");
1993 const Epetra_MultiVector &src,
1994 const Epetra_MultiVector &dst,
1999 Assert(src.Map().SameAs(mtrx.DomainMap()) ==
true,
2001 "Column map of matrix does not fit with vector map!"));
2002 Assert(dst.Map().SameAs(mtrx.RangeMap()) ==
true,
2003 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2007 Assert(src.Map().SameAs(mtrx.RangeMap()) ==
true,
2009 "Column map of matrix does not fit with vector map!"));
2010 Assert(dst.Map().SameAs(mtrx.DomainMap()) ==
true,
2011 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2020 const Epetra_MultiVector &src,
2021 const Epetra_MultiVector &dst,
2026 Assert(src.Map().SameAs(op.OperatorDomainMap()) ==
true,
2028 "Column map of operator does not fit with vector map!"));
2029 Assert(dst.Map().SameAs(op.OperatorRangeMap()) ==
true,
2031 "Row map of operator does not fit with vector map!"));
2035 Assert(src.Map().SameAs(op.OperatorRangeMap()) ==
true,
2037 "Column map of operator does not fit with vector map!"));
2038 Assert(dst.Map().SameAs(op.OperatorDomainMap()) ==
true,
2040 "Row map of operator does not fit with vector map!"));
2047 namespace LinearOperatorImplementation
2187 template <
typename Solver,
typename Preconditioner>
2189 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2191 Preconditioner>::value,
2212 template <
typename Solver,
typename Preconditioner>
2214 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2216 Preconditioner>::value),
2377 virtual const char *
2378 Label()
const override;
2387 virtual const Epetra_Comm &
2388 Comm()
const override;
2397 virtual const Epetra_Map &
2408 virtual const Epetra_Map &
2422 template <
typename EpetraOpType>
2424 const bool supports_inverse_operations,
2510 visit_present_row();
2514 inline AccessorBase::size_type
2515 AccessorBase::row()
const
2517 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2522 inline AccessorBase::size_type
2523 AccessorBase::column()
const
2526 return (*colnum_cache)[a_index];
2530 inline AccessorBase::size_type
2531 AccessorBase::index()
const
2538 inline Accessor<true>::Accessor(MatrixType * matrix,
2539 const size_type row,
2540 const size_type index)
2545 template <
bool Other>
2546 inline Accessor<true>::Accessor(
const Accessor<Other> &other)
2547 : AccessorBase(other)
2552 Accessor<true>::value()
const
2555 return (*value_cache)[a_index];
2559 inline Accessor<false>::Reference::Reference(
const Accessor<false> &acc)
2560 : accessor(const_cast<Accessor<false> &>(acc))
2564 inline Accessor<false>::Reference::operator
TrilinosScalar()
const
2566 return (*accessor.value_cache)[accessor.a_index];
2569 inline const Accessor<false>::Reference &
2570 Accessor<false>::Reference::operator=(
const TrilinosScalar n)
const
2572 (*accessor.value_cache)[accessor.a_index] = n;
2573 accessor.matrix->set(accessor.row(),
2580 inline const Accessor<false>::Reference &
2581 Accessor<false>::Reference::operator+=(
const TrilinosScalar n)
const
2583 (*accessor.value_cache)[accessor.a_index] += n;
2584 accessor.matrix->set(accessor.row(),
2591 inline const Accessor<false>::Reference &
2592 Accessor<false>::Reference::operator-=(
const TrilinosScalar n)
const
2594 (*accessor.value_cache)[accessor.a_index] -= n;
2595 accessor.matrix->set(accessor.row(),
2602 inline const Accessor<false>::Reference &
2603 Accessor<false>::Reference::operator*=(
const TrilinosScalar n)
const
2605 (*accessor.value_cache)[accessor.a_index] *= n;
2606 accessor.matrix->set(accessor.row(),
2613 inline const Accessor<false>::Reference &
2614 Accessor<false>::Reference::operator/=(
const TrilinosScalar n)
const
2616 (*accessor.value_cache)[accessor.a_index] /= n;
2617 accessor.matrix->set(accessor.row(),
2624 inline Accessor<false>::Accessor(MatrixType * matrix,
2625 const size_type row,
2626 const size_type index)
2631 inline Accessor<false>::Reference
2632 Accessor<false>::value()
const
2640 template <
bool Constness>
2641 inline Iterator<Constness>::Iterator(MatrixType * matrix,
2642 const size_type row,
2643 const size_type index)
2648 template <
bool Constness>
2649 template <
bool Other>
2650 inline Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2651 : accessor(other.accessor)
2655 template <
bool Constness>
2656 inline Iterator<Constness> &
2657 Iterator<Constness>::operator++()
2667 if (accessor.a_index >= accessor.colnum_cache->size())
2669 accessor.a_index = 0;
2672 while ((accessor.a_row < accessor.matrix->m()) &&
2673 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2674 (accessor.matrix->row_length(accessor.a_row) == 0)))
2677 accessor.visit_present_row();
2683 template <
bool Constness>
2684 inline Iterator<Constness>
2685 Iterator<Constness>::operator++(
int)
2687 const Iterator<Constness> old_state = *
this;
2694 template <
bool Constness>
2695 inline const Accessor<Constness> &
2696 Iterator<Constness>::operator*()
const
2703 template <
bool Constness>
2704 inline const Accessor<Constness> *
2705 Iterator<Constness>::operator->()
const
2712 template <
bool Constness>
2713 template <
bool OtherConstness>
2715 Iterator<Constness>::operator==(
const Iterator<OtherConstness> &other)
const
2717 return (accessor.a_row == other.accessor.a_row &&
2718 accessor.a_index == other.accessor.a_index);
2723 template <
bool Constness>
2724 template <
bool OtherConstness>
2726 Iterator<Constness>::operator!=(
const Iterator<OtherConstness> &other)
const
2728 return !(*
this == other);
2733 template <
bool Constness>
2734 template <
bool OtherConstness>
2736 Iterator<Constness>::operator<(
const Iterator<OtherConstness> &other)
const
2738 return (accessor.row() < other.accessor.row() ||
2739 (accessor.row() == other.accessor.row() &&
2740 accessor.index() < other.accessor.index()));
2744 template <
bool Constness>
2745 template <
bool OtherConstness>
2747 Iterator<Constness>::operator>(
const Iterator<OtherConstness> &other)
const
2749 return (other < *
this);
2767 return const_iterator(
this, m(), 0);
2776 if (in_local_range(r) && (row_length(r) > 0))
2777 return const_iterator(
this, r, 0);
2792 for (size_type i = r + 1; i < m(); ++i)
2793 if (in_local_range(i) && (row_length(i) > 0))
2794 return const_iterator(
this, i, 0);
2814 return iterator(
this, m(), 0);
2823 if (in_local_range(r) && (row_length(r) > 0))
2824 return iterator(
this, r, 0);
2839 for (size_type i = r + 1; i < m(); ++i)
2840 if (in_local_range(i) && (row_length(i) > 0))
2841 return iterator(
this, i, 0);
2854# ifndef DEAL_II_WITH_64BIT_INDICES
2859 end =
matrix->RowMap().MaxMyGID64() + 1;
2881 SparseMatrix::set<TrilinosScalar>(
const size_type row,
2882 const size_type n_cols,
2883 const size_type * col_indices,
2885 const bool elide_zero_values);
2889 template <
typename Number>
2892 const size_type n_cols,
2893 const size_type *col_indices,
2894 const Number * values,
2895 const bool elide_zero_values)
2897 std::vector<TrilinosScalar> trilinos_values(n_cols);
2898 std::copy(values, values + n_cols, trilinos_values.begin());
2900 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
2912 set(i, 1, &j, &value,
false);
2920 const bool elide_zero_values)
2926 for (size_type i = 0; i < indices.size(); ++i)
2951 if (last_action == Insert)
2953 const int ierr =
matrix->GlobalAssemble(*column_space_map,
2965 add(i, 1, &j, &value,
false);
2975# ifndef DEAL_II_WITH_64BIT_INDICES
2976 return matrix->NumGlobalRows();
2978 return matrix->NumGlobalRows64();
2999 return matrix->NumMyRows();
3004 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
3008# ifndef DEAL_II_WITH_64BIT_INDICES
3013 end =
matrix->RowMap().MaxMyGID64() + 1;
3021 inline std::uint64_t
3028 return static_cast<std::uint64_t
>(
matrix->NumGlobalNonzeros64());
3033 template <
typename SparsityPatternType>
3036 const SparsityPatternType &sparsity_pattern,
3038 const bool exchange_data)
3040 reinit(parallel_partitioning,
3041 parallel_partitioning,
3049 template <
typename number>
3052 const ::SparseMatrix<number> &sparse_matrix,
3054 const double drop_tolerance,
3055 const bool copy_values,
3056 const ::SparsityPattern * use_this_sparsity)
3060 reinit(parallel_partitioning,
3061 parallel_partitioning,
3070 inline const Epetra_CrsMatrix &
3073 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
3078 inline const Epetra_CrsGraph &
3119 namespace LinearOperatorImplementation
3121 template <
typename EpetraOpType>
3122 TrilinosPayload::TrilinosPayload(
3124 const bool supports_inverse_operations,
3125 const bool use_transpose,
3127 const IndexSet &locally_owned_domain_indices,
3128 const IndexSet &locally_owned_range_indices)
3129 : use_transpose(use_transpose)
3130 , communicator(mpi_communicator)
3132 locally_owned_domain_indices.make_trilinos_map(communicator.Comm()))
3134 locally_owned_range_indices.make_trilinos_map(communicator.Comm()))
3139 Assert(&tril_src != &tril_dst,
3146 const int ierr = op.Apply(tril_src, tril_dst);
3153 Assert(&tril_src != &tril_dst,
3158 !op.UseTranspose());
3160 op.SetUseTranspose(!op.UseTranspose());
3161 const int ierr = op.Apply(tril_src, tril_dst);
3163 op.SetUseTranspose(!op.UseTranspose());
3166 if (supports_inverse_operations)
3172 &tril_src != &tril_dst,
3177 !op.UseTranspose());
3179 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3187 &tril_src != &tril_dst,
3194 op.SetUseTranspose(!op.UseTranspose());
3195 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3197 op.SetUseTranspose(!op.UseTranspose());
3205 "Uninitialized TrilinosPayload::inv_vmult called. "
3206 "The operator does not support inverse operations."));
3212 "Uninitialized TrilinosPayload::inv_Tvmult called. "
3213 "The operator does not support inverse operations."));
3219 template <
typename Solver,
typename Preconditioner>
3221 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3223 Preconditioner>::value,
3227 const Preconditioner &preconditioner)
const
3229 const auto &payload = *
this;
3235 return_op.inv_vmult = [payload, &solver, &preconditioner](
3240 Assert(&tril_src != &tril_dst,
3245 !payload.UseTranspose());
3246 solver.solve(payload, tril_dst, tril_src, preconditioner);
3249 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3254 Assert(&tril_src != &tril_dst,
3259 payload.UseTranspose());
3262 solver.solve(payload, tril_dst, tril_src, preconditioner);
3268 if (return_op.UseTranspose() ==
true)
3269 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3274 template <
typename Solver,
typename Preconditioner>
3276 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3278 Preconditioner>::value),
3282 TrilinosPayload return_op(*
this);
3287 ExcMessage(
"Payload inv_vmult disabled because of "
3288 "incompatible solver/preconditioner choice."));
3294 ExcMessage(
"Payload inv_vmult disabled because of "
3295 "incompatible solver/preconditioner choice."));
3305 SparseMatrix::set<TrilinosScalar>(
const size_type row,
3306 const size_type n_cols,
3307 const size_type * col_indices,
3309 const bool elide_zero_values);
Epetra_Map make_trilinos_map(const MPI_Comm communicator=MPI_COMM_WORLD, const bool overlapping=false) const
AccessorBase(SparseMatrix *matrix, const size_type row, const size_type index)
std::shared_ptr< std::vector< TrilinosScalar > > value_cache
std::shared_ptr< std::vector< size_type > > colnum_cache
Reference(const Accessor< false > &accessor)
const Reference & operator-=(const TrilinosScalar n) const
const Reference & operator*=(const TrilinosScalar n) const
const Reference & operator+=(const TrilinosScalar n) const
const Reference & operator/=(const TrilinosScalar n) const
const Reference & operator=(const TrilinosScalar n) const
Accessor(MatrixType *matrix, const size_type row, const size_type index)
Accessor(const Accessor< Other > &a)
Accessor(MatrixType *matrix, const size_type row, const size_type index)
TrilinosScalar value() const
TrilinosScalar value() const
const Accessor< Constness > & operator*() const
const Accessor< Constness > * operator->() const
typename Accessor< Constness >::MatrixType MatrixType
bool operator==(const Iterator< OtherConstness > &) const
bool operator!=(const Iterator< OtherConstness > &) const
bool operator<(const Iterator< OtherConstness > &) const
Iterator(const Iterator< Other > &other)
::types::global_dof_index size_type
TrilinosScalar value_type
bool operator>(const Iterator< OtherConstness > &) const
Iterator(MatrixType *matrix, const size_type row, const size_type index)
Accessor< Constness > accessor
Iterator< Constness > & operator++()
Iterator< Constness > operator++(int)
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::unique_ptr< Epetra_Map > column_space_map
MPI_Comm get_mpi_communicator() const
SparseMatrix & operator*=(const TrilinosScalar factor)
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
TrilinosScalar l1_norm() const
void compress(VectorOperation::values operation)
std::unique_ptr< Epetra_FECrsMatrix > matrix
::types::global_dof_index size_type
TrilinosScalar linfty_norm() const
const Epetra_CrsMatrix & trilinos_matrix() const
size_type memory_consumption() const
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
Epetra_CombineMode last_action
IndexSet locally_owned_range_indices() const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
const_iterator begin() const
void reinit(const SparsityPatternType &sparsity_pattern)
std::enable_if_t< std::is_same< typename VectorType::value_type, TrilinosScalar >::value > Tvmult(VectorType &dst, const VectorType &src) const
void vmult_add(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const SparseMatrix &)=delete
const_iterator end(const size_type r) const
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
bool is_compressed() const
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() const
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
bool in_local_range(const size_type index) const
void copy_from(const SparseMatrix &source)
unsigned int row_length(const size_type row) const
TrilinosScalar frobenius_norm() const
std::uint64_t n_nonzero_elements() const
const_iterator begin(const size_type r) const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
SparseMatrix(const SparseMatrix &)=delete
iterator end(const size_type r)
TrilinosScalar diag_element(const size_type i) const
void add(const size_type i, const size_type j, const TrilinosScalar value)
unsigned int local_size() const
iterator begin(const size_type r)
TrilinosScalar operator()(const size_type i, const size_type j) const
std::enable_if_t< !std::is_same< SparsityPatternType, ::SparseMatrix< double > >::value > 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)
void reinit(const IndexSet ¶llel_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)
virtual ~SparseMatrix() override=default
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)
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
const_iterator end() const
void set(const std::vector< size_type > &indices, const FullMatrix< TrilinosScalar > &full_matrix, const bool elide_zero_values=false)
std::pair< size_type, size_type > local_range() const
std::enable_if_t< std::is_same< typename VectorType::value_type, TrilinosScalar >::value > vmult(VectorType &dst, const VectorType &src) const
void reinit(const IndexSet ¶llel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm communicator=MPI_COMM_WORLD, const bool exchange_data=false)
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
TrilinosScalar value_type
Epetra_MpiComm communicator
std::enable_if_t< !(std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value), TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
virtual int SetUseTranspose(bool UseTranspose) override
virtual bool UseTranspose() const override
virtual const Epetra_Map & OperatorDomainMap() const override
IndexSet locally_owned_range_indices() const
TrilinosPayload transpose_payload() const
Epetra_MultiVector VectorType
MPI_Comm get_mpi_communicator() const
virtual ~TrilinosPayload() override=default
virtual int ApplyInverse(const VectorType &Y, VectorType &X) const override
std::function< void(VectorType &, const VectorType &)> Tvmult
virtual int Apply(const VectorType &X, VectorType &Y) const override
std::function< void(VectorType &, const VectorType &)> vmult
virtual const Epetra_Map & OperatorRangeMap() const override
std::function< void(VectorType &, const VectorType &)> inv_vmult
virtual const char * Label() const override
IndexSet locally_owned_domain_indices() const
virtual const Epetra_Comm & Comm() const override
virtual bool HasNormInf() const override
TrilinosPayload identity_payload() const
std::enable_if_t< std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value, TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
TrilinosPayload(EpetraOpType &op, const bool supports_inverse_operations, const bool use_transpose, const MPI_Comm mpi_communicator, const IndexSet &locally_owned_domain_indices, const IndexSet &locally_owned_range_indices)
virtual double NormInf() const override
TrilinosPayload null_payload() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
__global__ void set(Number *val, const Number s, const size_type N)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
#define DeclException0(Exception0)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcIteratorPastEnd()
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
#define AssertIsFinite(number)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
static ::ExceptionBase & ExcBeyondEndOfMatrix()
#define DeclException2(Exception2, type1, type2, outsequence)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcSourceEqualsDestination()
#define DeclException3(Exception3, type1, type2, type3, outsequence)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcNotQuadratic()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcMatrixNotCompressed()
static ::ExceptionBase & ExcTrilinosError(int arg1)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define AssertThrow(cond, exc)
@ matrix
Contents is actually a matrix.
types::global_dof_index size_type
TrilinosPayload operator+(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
TrilinosPayload operator*(const TrilinosPayload &first_op, const TrilinosPayload &second_op)
VectorType::value_type * end(VectorType &V)
VectorType::value_type * begin(VectorType &V)
void check_vector_map_equality(const Epetra_CrsMatrix &mtrx, const Epetra_MultiVector &src, const Epetra_MultiVector &dst, const bool transpose)
TrilinosWrappers::types::int64_type n_global_elements(const Epetra_BlockMap &map)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
unsigned int global_dof_index
static const bool zero_addition_can_be_elided
typename ::TrilinosWrappers::SparseMatrixIterators::Iterator< Constness >::value_type value_type
forward_iterator_tag iterator_category
typename ::TrilinosWrappers::SparseMatrixIterators::Iterator< Constness >::difference_type difference_type