16#ifndef dealii_trilinos_sparse_matrix_h
17# define dealii_trilinos_sparse_matrix_h
22# ifdef DEAL_II_WITH_TRILINOS
36# include <Epetra_Comm.h>
37# include <Epetra_CrsGraph.h>
38# include <Epetra_Export.h>
39# include <Epetra_FECrsMatrix.h>
40# include <Epetra_Map.h>
41# include <Epetra_MpiComm.h>
42# include <Epetra_MultiVector.h>
43# 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>
447 <<
"Attempt to access element " << arg2 <<
" of row "
448 << arg1 <<
" which doesn't have that many elements.");
456 template <
bool Other>
467 template <
bool Constness>
468 struct iterator_traits<
473 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
474 Constness>::value_type;
476 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
477 Constness>::difference_type;
559 <<
"You tried to access row " << arg1
560 <<
" of a non-contiguous locally owned row set."
561 <<
" The row " << arg1
562 <<
" is not stored locally and can't be accessed.");
613 const unsigned int n_max_entries_per_row);
624 const std::vector<unsigned int> &n_entries_per_row);
668 template <
typename SparsityPatternType>
670 reinit(
const SparsityPatternType &sparsity_pattern);
718 template <
typename number>
720 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
721 const double drop_tolerance = 1e-13,
722 const bool copy_values =
true,
723 const ::SparsityPattern * use_this_sparsity =
nullptr);
731 reinit(
const Epetra_CrsMatrix &input_matrix,
const bool copy_values =
true);
751 const MPI_Comm & communicator = MPI_COMM_WORLD,
752 const unsigned int n_max_entries_per_row = 0);
763 const std::vector<unsigned int> &n_entries_per_row);
780 const IndexSet &col_parallel_partitioning,
781 const MPI_Comm &communicator = MPI_COMM_WORLD,
782 const size_type n_max_entries_per_row = 0);
799 const IndexSet & col_parallel_partitioning,
801 const std::vector<unsigned int> &n_entries_per_row);
823 template <
typename SparsityPatternType>
826 const SparsityPatternType &sparsity_pattern,
827 const MPI_Comm & communicator = MPI_COMM_WORLD,
828 const bool exchange_data =
false);
842 template <
typename SparsityPatternType>
843 typename std::enable_if<
844 !std::is_same<SparsityPatternType,
847 const IndexSet & col_parallel_partitioning,
848 const SparsityPatternType &sparsity_pattern,
849 const MPI_Comm & communicator = MPI_COMM_WORLD,
850 const bool exchange_data =
false);
868 template <
typename number>
871 const ::SparseMatrix<number> &dealii_sparse_matrix,
872 const MPI_Comm & communicator = MPI_COMM_WORLD,
873 const double drop_tolerance = 1e-13,
874 const bool copy_values =
true,
875 const ::SparsityPattern * use_this_sparsity =
nullptr);
890 template <
typename number>
893 const IndexSet & col_parallel_partitioning,
894 const ::SparseMatrix<number> &dealii_sparse_matrix,
895 const MPI_Comm & communicator = MPI_COMM_WORLD,
896 const double drop_tolerance = 1e-13,
897 const bool copy_values =
true,
898 const ::SparsityPattern * use_this_sparsity =
nullptr);
936 std::pair<size_type, size_type>
1097 set(
const std::vector<size_type> & indices,
1099 const bool elide_zero_values =
false);
1107 set(
const std::vector<size_type> & row_indices,
1108 const std::vector<size_type> & col_indices,
1110 const bool elide_zero_values =
false);
1141 const std::vector<size_type> & col_indices,
1142 const std::vector<TrilinosScalar> &values,
1143 const bool elide_zero_values =
false);
1172 template <
typename Number>
1177 const Number * values,
1178 const bool elide_zero_values =
false);
1211 add(
const std::vector<size_type> & indices,
1213 const bool elide_zero_values =
true);
1221 add(
const std::vector<size_type> & row_indices,
1222 const std::vector<size_type> & col_indices,
1224 const bool elide_zero_values =
true);
1241 const std::vector<size_type> & col_indices,
1242 const std::vector<TrilinosScalar> &values,
1243 const bool elide_zero_values =
true);
1263 const bool elide_zero_values =
true,
1264 const bool col_indices_are_sorted =
false);
1344 clear_rows(
const std::vector<size_type> &rows,
1440 template <
typename VectorType>
1441 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1443 vmult(VectorType &dst,
const VectorType &src)
const;
1451 template <
typename VectorType>
1452 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1454 vmult(VectorType &dst,
const VectorType &src)
const;
1470 template <
typename VectorType>
1471 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1473 Tvmult(VectorType &dst,
const VectorType &src)
const;
1481 template <
typename VectorType>
1482 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1484 Tvmult(VectorType &dst,
const VectorType &src)
const;
1495 template <
typename VectorType>
1497 vmult_add(VectorType &dst,
const VectorType &src)
const;
1509 template <
typename VectorType>
1511 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
1665 const Epetra_CrsMatrix &
1672 const Epetra_CrsGraph &
1821 print(std::ostream &out,
1822 const bool write_extended_trilinos_info =
false)
const;
1835 <<
"An error with error number " << arg1
1836 <<
" occurred while calling a Trilinos function. "
1838 "For historical reasons, many Trilinos functions express "
1839 "errors by returning specific integer values to indicate "
1840 "certain errors. Unfortunately, different Trilinos functions "
1841 "often use the same integer values for different kinds of "
1842 "errors, and in most cases it is also not documented what "
1843 "each error code actually means. As a consequence, it is often "
1844 "difficult to say what a particular error (in this case, "
1845 "the error with integer code '"
1847 <<
"') represents and how one should fix a code to avoid it. "
1848 "The best one can often do is to look up the call stack to "
1849 "see which deal.II function generated the error, and which "
1850 "Trilinos function the error code had originated from; "
1851 "then look up the Trilinos source code of that function (for "
1852 "example on github) to see what code path set that error "
1853 "code. Short of going through all of that, the only other "
1854 "option is to guess the cause of the error from "
1855 "the context in which the error appeared.");
1864 <<
"The entry with index <" << arg1 <<
',' << arg2
1865 <<
"> does not exist.");
1871 "You are attempting an operation on two matrices that "
1872 "are the same object, but the operation requires that the "
1873 "two objects are in fact different.");
1888 <<
"You tried to access element (" << arg1 <<
'/' << arg2
1890 <<
" of a distributed matrix, but only rows in range ["
1891 << arg3 <<
',' << arg4
1892 <<
"] are stored locally and can be accessed.");
1900 <<
"You tried to access element (" << arg1 <<
'/' << arg2
1901 <<
')' <<
" of a sparse matrix, but it appears to not"
1902 <<
" exist in the Trilinos sparsity pattern.");
1990 const Epetra_MultiVector &src,
1991 const Epetra_MultiVector &dst,
1996 Assert(src.Map().SameAs(mtrx.DomainMap()) ==
true,
1998 "Column map of matrix does not fit with vector map!"));
1999 Assert(dst.Map().SameAs(mtrx.RangeMap()) ==
true,
2000 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2004 Assert(src.Map().SameAs(mtrx.RangeMap()) ==
true,
2006 "Column map of matrix does not fit with vector map!"));
2007 Assert(dst.Map().SameAs(mtrx.DomainMap()) ==
true,
2008 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2017 const Epetra_MultiVector &src,
2018 const Epetra_MultiVector &dst,
2023 Assert(src.Map().SameAs(op.OperatorDomainMap()) ==
true,
2025 "Column map of operator does not fit with vector map!"));
2026 Assert(dst.Map().SameAs(op.OperatorRangeMap()) ==
true,
2028 "Row map of operator does not fit with vector map!"));
2032 Assert(src.Map().SameAs(op.OperatorRangeMap()) ==
true,
2034 "Column map of operator does not fit with vector map!"));
2035 Assert(dst.Map().SameAs(op.OperatorDomainMap()) ==
true,
2037 "Row map of operator does not fit with vector map!"));
2044 namespace LinearOperatorImplementation
2184 template <
typename Solver,
typename Preconditioner>
2185 typename std::enable_if<
2186 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2188 Preconditioner>::value,
2209 template <
typename Solver,
typename Preconditioner>
2210 typename std::enable_if<
2211 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2213 Preconditioner>::value),
2374 virtual const char *
2375 Label()
const override;
2384 virtual const Epetra_Comm &
2385 Comm()
const override;
2394 virtual const Epetra_Map &
2405 virtual const Epetra_Map &
2419 template <
typename EpetraOpType>
2421 const bool supports_inverse_operations,
2507 visit_present_row();
2511 inline AccessorBase::size_type
2512 AccessorBase::row()
const
2514 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2519 inline AccessorBase::size_type
2520 AccessorBase::column()
const
2523 return (*colnum_cache)[a_index];
2527 inline AccessorBase::size_type
2528 AccessorBase::index()
const
2535 inline Accessor<true>::Accessor(MatrixType * matrix,
2536 const size_type row,
2537 const size_type index)
2542 template <
bool Other>
2543 inline Accessor<true>::Accessor(
const Accessor<Other> &other)
2544 : AccessorBase(other)
2549 Accessor<true>::value()
const
2552 return (*value_cache)[a_index];
2556 inline Accessor<false>::Reference::Reference(
const Accessor<false> &acc)
2557 : accessor(const_cast<Accessor<false> &>(acc))
2561 inline Accessor<false>::Reference::operator
TrilinosScalar()
const
2563 return (*accessor.value_cache)[accessor.a_index];
2566 inline const Accessor<false>::Reference &
2567 Accessor<false>::Reference::operator=(
const TrilinosScalar n)
const
2569 (*accessor.value_cache)[accessor.a_index] = n;
2570 accessor.matrix->set(accessor.row(),
2577 inline const Accessor<false>::Reference &
2578 Accessor<false>::Reference::operator+=(
const TrilinosScalar n)
const
2580 (*accessor.value_cache)[accessor.a_index] += n;
2581 accessor.matrix->set(accessor.row(),
2588 inline const Accessor<false>::Reference &
2589 Accessor<false>::Reference::operator-=(
const TrilinosScalar n)
const
2591 (*accessor.value_cache)[accessor.a_index] -= n;
2592 accessor.matrix->set(accessor.row(),
2599 inline const Accessor<false>::Reference &
2600 Accessor<false>::Reference::operator*=(
const TrilinosScalar n)
const
2602 (*accessor.value_cache)[accessor.a_index] *= n;
2603 accessor.matrix->set(accessor.row(),
2610 inline const Accessor<false>::Reference &
2611 Accessor<false>::Reference::operator/=(
const TrilinosScalar n)
const
2613 (*accessor.value_cache)[accessor.a_index] /= n;
2614 accessor.matrix->set(accessor.row(),
2621 inline Accessor<false>::Accessor(MatrixType * matrix,
2622 const size_type row,
2623 const size_type index)
2628 inline Accessor<false>::Reference
2629 Accessor<false>::value()
const
2637 template <
bool Constness>
2638 inline Iterator<Constness>::Iterator(MatrixType * matrix,
2639 const size_type row,
2640 const size_type index)
2645 template <
bool Constness>
2646 template <
bool Other>
2647 inline Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2648 : accessor(other.accessor)
2652 template <
bool Constness>
2653 inline Iterator<Constness> &
2654 Iterator<Constness>::operator++()
2664 if (accessor.a_index >= accessor.colnum_cache->size())
2666 accessor.a_index = 0;
2669 while ((accessor.a_row < accessor.matrix->m()) &&
2670 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2671 (accessor.matrix->row_length(accessor.a_row) == 0)))
2674 accessor.visit_present_row();
2680 template <
bool Constness>
2681 inline Iterator<Constness>
2682 Iterator<Constness>::operator++(
int)
2684 const Iterator<Constness> old_state = *
this;
2691 template <
bool Constness>
2692 inline const Accessor<Constness> &
2693 Iterator<Constness>::operator*()
const
2700 template <
bool Constness>
2701 inline const Accessor<Constness> *
2702 Iterator<Constness>::operator->()
const
2709 template <
bool Constness>
2711 Iterator<Constness>::operator==(
const Iterator<Constness> &other)
const
2713 return (accessor.a_row == other.accessor.a_row &&
2714 accessor.a_index == other.accessor.a_index);
2719 template <
bool Constness>
2721 Iterator<Constness>::operator!=(
const Iterator<Constness> &other)
const
2723 return !(*
this == other);
2728 template <
bool Constness>
2730 Iterator<Constness>::operator<(
const Iterator<Constness> &other)
const
2732 return (accessor.row() < other.accessor.row() ||
2733 (accessor.row() == other.accessor.row() &&
2734 accessor.index() < other.accessor.index()));
2738 template <
bool Constness>
2740 Iterator<Constness>::operator>(
const Iterator<Constness> &other)
const
2742 return (other < *
this);
2760 return const_iterator(
this, m(), 0);
2769 if (in_local_range(r) && (row_length(r) > 0))
2770 return const_iterator(
this, r, 0);
2785 for (size_type i = r + 1; i < m(); ++i)
2786 if (in_local_range(i) && (row_length(i) > 0))
2787 return const_iterator(
this, i, 0);
2807 return iterator(
this, m(), 0);
2816 if (in_local_range(r) && (row_length(r) > 0))
2817 return iterator(
this, r, 0);
2832 for (size_type i = r + 1; i < m(); ++i)
2833 if (in_local_range(i) && (row_length(i) > 0))
2834 return iterator(
this, i, 0);
2847# ifndef DEAL_II_WITH_64BIT_INDICES
2852 end =
matrix->RowMap().MaxMyGID64() + 1;
2874 SparseMatrix::set<TrilinosScalar>(
const size_type row,
2875 const size_type n_cols,
2876 const size_type * col_indices,
2878 const bool elide_zero_values);
2882 template <
typename Number>
2885 const size_type n_cols,
2886 const size_type *col_indices,
2887 const Number * values,
2888 const bool elide_zero_values)
2890 std::vector<TrilinosScalar> trilinos_values(n_cols);
2891 std::copy(values, values + n_cols, trilinos_values.begin());
2893 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
2905 set(i, 1, &j, &value,
false);
2913 const bool elide_zero_values)
2919 for (size_type i = 0; i < indices.size(); ++i)
2944 if (last_action == Insert)
2946 const int ierr =
matrix->GlobalAssemble(*column_space_map,
2958 add(i, 1, &j, &value,
false);
2968# ifndef DEAL_II_WITH_64BIT_INDICES
2969 return matrix->NumGlobalRows();
2971 return matrix->NumGlobalRows64();
2992 return matrix->NumMyRows();
2997 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
3001# ifndef DEAL_II_WITH_64BIT_INDICES
3006 end =
matrix->RowMap().MaxMyGID64() + 1;
3014 inline std::uint64_t
3021 return static_cast<std::uint64_t
>(
matrix->NumGlobalNonzeros64());
3026 template <
typename SparsityPatternType>
3029 const SparsityPatternType &sparsity_pattern,
3031 const bool exchange_data)
3033 reinit(parallel_partitioning,
3034 parallel_partitioning,
3042 template <
typename number>
3045 const ::SparseMatrix<number> &sparse_matrix,
3047 const double drop_tolerance,
3048 const bool copy_values,
3049 const ::SparsityPattern * use_this_sparsity)
3053 reinit(parallel_partitioning,
3054 parallel_partitioning,
3063 inline const Epetra_CrsMatrix &
3066 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
3071 inline const Epetra_CrsGraph &
3112 namespace LinearOperatorImplementation
3114 template <
typename EpetraOpType>
3115 TrilinosPayload::TrilinosPayload(
3117 const bool supports_inverse_operations,
3118 const bool use_transpose,
3120 const IndexSet &locally_owned_domain_indices,
3121 const IndexSet &locally_owned_range_indices)
3122 : use_transpose(use_transpose)
3123 , communicator(mpi_communicator)
3125 locally_owned_domain_indices.make_trilinos_map(communicator.Comm()))
3127 locally_owned_range_indices.make_trilinos_map(communicator.Comm()))
3132 Assert(&tril_src != &tril_dst,
3139 const int ierr = op.Apply(tril_src, tril_dst);
3146 Assert(&tril_src != &tril_dst,
3151 !op.UseTranspose());
3153 op.SetUseTranspose(!op.UseTranspose());
3154 const int ierr = op.Apply(tril_src, tril_dst);
3156 op.SetUseTranspose(!op.UseTranspose());
3159 if (supports_inverse_operations)
3165 &tril_src != &tril_dst,
3170 !op.UseTranspose());
3172 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3180 &tril_src != &tril_dst,
3187 op.SetUseTranspose(!op.UseTranspose());
3188 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3190 op.SetUseTranspose(!op.UseTranspose());
3198 "Uninitialized TrilinosPayload::inv_vmult called. "
3199 "The operator does not support inverse operations."));
3205 "Uninitialized TrilinosPayload::inv_Tvmult called. "
3206 "The operator does not support inverse operations."));
3212 template <
typename Solver,
typename Preconditioner>
3213 typename std::enable_if<
3214 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3216 Preconditioner>::value,
3217 TrilinosPayload>::type
3220 const Preconditioner &preconditioner)
const
3222 const auto &payload = *
this;
3228 return_op.inv_vmult = [payload, &solver, &preconditioner](
3233 Assert(&tril_src != &tril_dst,
3238 !payload.UseTranspose());
3239 solver.solve(payload, tril_dst, tril_src, preconditioner);
3242 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3247 Assert(&tril_src != &tril_dst,
3252 payload.UseTranspose());
3255 solver.solve(payload, tril_dst, tril_src, preconditioner);
3261 if (return_op.UseTranspose() ==
true)
3262 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3267 template <
typename Solver,
typename Preconditioner>
3268 typename std::enable_if<
3269 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3271 Preconditioner>::value),
3272 TrilinosPayload>::type
3275 TrilinosPayload return_op(*
this);
3280 ExcMessage(
"Payload inv_vmult disabled because of "
3281 "incompatible solver/preconditioner choice."));
3287 ExcMessage(
"Payload inv_vmult disabled because of "
3288 "incompatible solver/preconditioner choice."));
3298 SparseMatrix::set<TrilinosScalar>(
const size_type row,
3299 const size_type n_cols,
3300 const size_type * col_indices,
3302 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
bool operator>(const Iterator< Constness > &) const
typename Accessor< Constness >::MatrixType MatrixType
Iterator(const Iterator< Other > &other)
::types::global_dof_index size_type
bool operator==(const Iterator< Constness > &) const
TrilinosScalar value_type
bool operator!=(const Iterator< Constness > &) const
Iterator(MatrixType *matrix, const size_type row, const size_type index)
Accessor< Constness > accessor
Iterator< Constness > & operator++()
bool operator<(const Iterator< Constness > &) const
Iterator< Constness > operator++(int)
void set(const size_type i, const size_type j, const TrilinosScalar value)
std::enable_if<!std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
std::unique_ptr< Epetra_Map > column_space_map
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
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
std::unique_ptr< Epetra_FECrsMatrix > matrix
std::enable_if<!std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
::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
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)
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)
void vmult_add(VectorType &dst, const VectorType &src) const
SparseMatrix & operator=(const SparseMatrix &)=delete
const_iterator end(const size_type r) const
void compress(::VectorOperation::values operation)
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)
void reinit(const IndexSet ¶llel_partitioning, const SparsityPatternType &sparsity_pattern, const MPI_Comm &communicator=MPI_COMM_WORLD, const bool exchange_data=false)
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)
std::enable_if< std::is_same< typenameVectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
TrilinosScalar operator()(const size_type i, const size_type j) const
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::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
TrilinosScalar value_type
Epetra_MpiComm communicator
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
std::enable_if< std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value, TrilinosPayload >::type 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 ~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
std::enable_if<!(std::is_base_of< TrilinosWrappers::SolverBase, Solver >::value &&std::is_base_of< TrilinosWrappers::PreconditionBase, Preconditioner >::value), TrilinosPayload >::type inverse_payload(Solver &, const Preconditioner &) const
IndexSet locally_owned_domain_indices() const
virtual const Epetra_Comm & Comm() const override
virtual bool HasNormInf() const override
TrilinosPayload identity_payload() const
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