15#ifndef dealii_trilinos_sparse_matrix_h
16# define dealii_trilinos_sparse_matrix_h
21# 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>
50# include <type_traits>
57template <
typename MatrixType>
60template <
typename number>
72 template <
bool Constness>
97 <<
"You tried to access row " << arg1
98 <<
" of a distributed sparsity pattern, "
99 <<
" but only rows " << arg2 <<
" through " << arg3
100 <<
" are stored locally and can be accessed.");
198 template <
bool Constess>
239 template <
bool Other>
352 template <
bool Constness>
388 template <
bool Other>
419 template <
bool OtherConstness>
426 template <
bool OtherConstness>
435 template <
bool OtherConstness>
442 template <
bool OtherConstness>
452 <<
"Attempt to access element " << arg2 <<
" of row "
453 << arg1 <<
" which doesn't have that many elements.");
461 template <
bool Other>
472 template <
bool Constness>
478 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
481 typename ::TrilinosWrappers::SparseMatrixIterators::Iterator<
564 <<
"You tried to access row " << arg1
565 <<
" of a non-contiguous locally owned row set."
566 <<
" The row " << arg1
567 <<
" is not stored locally and can't be accessed.");
618 const unsigned int n_max_entries_per_row);
629 const std::vector<unsigned int> &n_entries_per_row);
673 template <
typename SparsityPatternType>
675 reinit(
const SparsityPatternType &sparsity_pattern);
723 template <
typename number>
725 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
726 const double drop_tolerance = 1e-13,
727 const bool copy_values =
true,
728 const ::SparsityPattern *use_this_sparsity =
nullptr);
736 reinit(
const Epetra_CrsMatrix &input_matrix,
const bool copy_values =
true);
756 const MPI_Comm communicator = MPI_COMM_WORLD,
757 const unsigned int n_max_entries_per_row = 0);
768 const std::vector<unsigned int> &n_entries_per_row);
785 const IndexSet &col_parallel_partitioning,
786 const MPI_Comm communicator = MPI_COMM_WORLD,
787 const size_type n_max_entries_per_row = 0);
804 const IndexSet &col_parallel_partitioning,
806 const std::vector<unsigned int> &n_entries_per_row);
828 template <
typename SparsityPatternType>
831 const SparsityPatternType &sparsity_pattern,
832 const MPI_Comm communicator = MPI_COMM_WORLD,
833 const bool exchange_data =
false);
847 template <
typename SparsityPatternType>
849 !std::is_same_v<SparsityPatternType, ::SparseMatrix<double>>>
851 const IndexSet &col_parallel_partitioning,
852 const SparsityPatternType &sparsity_pattern,
853 const MPI_Comm communicator = MPI_COMM_WORLD,
854 const bool exchange_data =
false);
872 template <
typename number>
875 const ::SparseMatrix<number> &dealii_sparse_matrix,
876 const MPI_Comm communicator = MPI_COMM_WORLD,
877 const double drop_tolerance = 1e-13,
878 const bool copy_values =
true,
879 const ::SparsityPattern *use_this_sparsity =
nullptr);
894 template <
typename number>
897 const IndexSet &col_parallel_partitioning,
898 const ::SparseMatrix<number> &dealii_sparse_matrix,
899 const MPI_Comm communicator = MPI_COMM_WORLD,
900 const double drop_tolerance = 1e-13,
901 const bool copy_values =
true,
902 const ::SparsityPattern *use_this_sparsity =
nullptr);
940 std::pair<size_type, size_type>
1101 set(
const std::vector<size_type> &indices,
1103 const bool elide_zero_values =
false);
1111 set(
const std::vector<size_type> &row_indices,
1112 const std::vector<size_type> &col_indices,
1114 const bool elide_zero_values =
false);
1145 const std::vector<size_type> &col_indices,
1146 const std::vector<TrilinosScalar> &values,
1147 const bool elide_zero_values =
false);
1176 template <
typename Number>
1181 const Number *values,
1182 const bool elide_zero_values =
false);
1215 add(
const std::vector<size_type> &indices,
1217 const bool elide_zero_values =
true);
1225 add(
const std::vector<size_type> &row_indices,
1226 const std::vector<size_type> &col_indices,
1228 const bool elide_zero_values =
true);
1245 const std::vector<size_type> &col_indices,
1246 const std::vector<TrilinosScalar> &values,
1247 const bool elide_zero_values =
true);
1267 const bool elide_zero_values =
true,
1268 const bool col_indices_are_sorted =
false);
1444 template <
typename VectorType>
1446 std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1447 vmult(VectorType &dst,
const VectorType &src)
const;
1455 template <
typename VectorType>
1457 !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1458 vmult(VectorType &dst,
const VectorType &src)
const;
1474 template <
typename VectorType>
1476 std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1477 Tvmult(VectorType &dst,
const VectorType &src)
const;
1485 template <
typename VectorType>
1487 !std::is_same_v<typename VectorType::value_type, TrilinosScalar>>
1488 Tvmult(VectorType &dst,
const VectorType &src)
const;
1499 template <
typename VectorType>
1501 vmult_add(VectorType &dst,
const VectorType &src)
const;
1513 template <
typename VectorType>
1515 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
1579 template <
typename VectorType>
1581 residual(VectorType &dst,
const VectorType &x,
const VectorType &b)
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.");
1995 const Epetra_MultiVector &src,
1996 const Epetra_MultiVector &dst,
2001 Assert(src.Map().SameAs(mtrx.DomainMap()) ==
true,
2003 "Column map of matrix does not fit with vector map!"));
2004 Assert(dst.Map().SameAs(mtrx.RangeMap()) ==
true,
2005 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2009 Assert(src.Map().SameAs(mtrx.RangeMap()) ==
true,
2011 "Column map of matrix does not fit with vector map!"));
2012 Assert(dst.Map().SameAs(mtrx.DomainMap()) ==
true,
2013 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2019 const Epetra_MultiVector &src,
2020 const Epetra_MultiVector &dst,
2025 Assert(src.Map().SameAs(op.OperatorDomainMap()) ==
true,
2027 "Column map of operator does not fit with vector map!"));
2028 Assert(dst.Map().SameAs(op.OperatorRangeMap()) ==
true,
2030 "Row map of operator does not fit with vector map!"));
2034 Assert(src.Map().SameAs(op.OperatorRangeMap()) ==
true,
2036 "Column map of operator does not fit with vector map!"));
2037 Assert(dst.Map().SameAs(op.OperatorDomainMap()) ==
true,
2039 "Row map of operator does not fit with vector map!"));
2189 template <
typename Solver,
typename Preconditioner>
2191 std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
2214 template <
typename Solver,
typename Preconditioner>
2216 !(std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
2379 virtual const char *
2380 Label()
const override;
2389 virtual const Epetra_Comm &
2390 Comm()
const override;
2399 virtual const Epetra_Map &
2410 virtual const Epetra_Map &
2424 template <
typename EpetraOpType>
2426 const bool supports_inverse_operations,
2512 visit_present_row();
2516 inline AccessorBase::size_type
2517 AccessorBase::row()
const
2519 Assert(a_row < matrix->m(), ExcBeyondEndOfMatrix());
2524 inline AccessorBase::size_type
2525 AccessorBase::column()
const
2528 return (*colnum_cache)[a_index];
2532 inline AccessorBase::size_type
2533 AccessorBase::index()
const
2540 inline Accessor<true>::Accessor(MatrixType *matrix,
2547 template <
bool Other>
2548 inline Accessor<true>::Accessor(
const Accessor<Other> &other)
2549 : AccessorBase(other)
2554 Accessor<true>::value()
const
2557 return (*value_cache)[a_index];
2561 inline Accessor<false>::Reference::Reference(
const Accessor<false> &acc)
2562 : accessor(const_cast<Accessor<false> &>(acc))
2566 inline Accessor<false>::Reference::operator
TrilinosScalar()
const
2568 return (*accessor.value_cache)[accessor.a_index];
2571 inline const Accessor<false>::Reference &
2572 Accessor<false>::Reference::operator=(
const TrilinosScalar n)
const
2574 (*accessor.value_cache)[accessor.a_index] = n;
2575 accessor.matrix->set(accessor.row(),
2582 inline const Accessor<false>::Reference &
2583 Accessor<false>::Reference::operator+=(
const TrilinosScalar n)
const
2585 (*accessor.value_cache)[accessor.a_index] += n;
2586 accessor.matrix->set(accessor.row(),
2593 inline const Accessor<false>::Reference &
2594 Accessor<false>::Reference::operator-=(
const TrilinosScalar n)
const
2596 (*accessor.value_cache)[accessor.a_index] -= n;
2597 accessor.matrix->set(accessor.row(),
2604 inline const Accessor<false>::Reference &
2605 Accessor<false>::Reference::operator*=(
const TrilinosScalar n)
const
2607 (*accessor.value_cache)[accessor.a_index] *= n;
2608 accessor.matrix->set(accessor.row(),
2615 inline const Accessor<false>::Reference &
2616 Accessor<false>::Reference::operator/=(
const TrilinosScalar n)
const
2618 (*accessor.value_cache)[accessor.a_index] /= n;
2619 accessor.matrix->set(accessor.row(),
2626 inline Accessor<false>::Accessor(MatrixType *matrix,
2633 inline Accessor<false>::Reference
2634 Accessor<false>::value()
const
2642 template <
bool Constness>
2643 inline Iterator<Constness>::Iterator(MatrixType *matrix,
2650 template <
bool Constness>
2651 template <
bool Other>
2652 inline Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2653 : accessor(other.accessor)
2657 template <
bool Constness>
2658 inline Iterator<Constness> &
2659 Iterator<Constness>::operator++()
2669 if (accessor.a_index >= accessor.colnum_cache->size())
2671 accessor.a_index = 0;
2674 while ((accessor.a_row < accessor.matrix->m()) &&
2675 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2676 (accessor.matrix->row_length(accessor.a_row) == 0)))
2679 accessor.visit_present_row();
2685 template <
bool Constness>
2686 inline Iterator<Constness>
2687 Iterator<Constness>::operator++(
int)
2689 const Iterator<Constness> old_state = *
this;
2696 template <
bool Constness>
2697 inline const Accessor<Constness> &
2698 Iterator<Constness>::operator*()
const
2705 template <
bool Constness>
2706 inline const Accessor<Constness> *
2707 Iterator<Constness>::operator->()
const
2714 template <
bool Constness>
2715 template <
bool OtherConstness>
2717 Iterator<Constness>::operator==(
const Iterator<OtherConstness> &other)
const
2719 return (accessor.a_row == other.accessor.a_row &&
2720 accessor.a_index == other.accessor.a_index);
2725 template <
bool Constness>
2726 template <
bool OtherConstness>
2728 Iterator<Constness>::operator!=(
const Iterator<OtherConstness> &other)
const
2730 return !(*
this == other);
2735 template <
bool Constness>
2736 template <
bool OtherConstness>
2738 Iterator<Constness>::operator<(
const Iterator<OtherConstness> &other)
const
2740 return (accessor.row() < other.accessor.row() ||
2741 (accessor.row() == other.accessor.row() &&
2742 accessor.index() < other.accessor.index()));
2746 template <
bool Constness>
2747 template <
bool OtherConstness>
2749 Iterator<Constness>::operator>(
const Iterator<OtherConstness> &other)
const
2751 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);
2791 for (size_type i = r + 1; i < m(); ++i)
2792 if (in_local_range(i) && (row_length(i) > 0))
2793 return const_iterator(
this, i, 0);
2811 return iterator(
this, m(), 0);
2819 if (in_local_range(r) && (row_length(r) > 0))
2820 return iterator(
this, r, 0);
2834 for (size_type i = r + 1; i < m(); ++i)
2835 if (in_local_range(i) && (row_length(i) > 0))
2836 return iterator(
this, i, 0);
2848# ifndef DEAL_II_WITH_64BIT_INDICES
2853 end =
matrix->RowMap().MaxMyGID64() + 1;
2874 const size_type n_cols,
2875 const size_type *col_indices,
2877 const bool elide_zero_values);
2881 template <
typename Number>
2883 const size_type n_cols,
2884 const size_type *col_indices,
2885 const Number *values,
2886 const bool elide_zero_values)
2888 std::vector<TrilinosScalar> trilinos_values(n_cols);
2889 std::copy(values, values + n_cols, trilinos_values.begin());
2891 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
2902 set(i, 1, &j, &
value,
false);
2908 const FullMatrix<TrilinosScalar> &values,
2909 const bool elide_zero_values)
2915 for (size_type i = 0; i < indices.size(); ++i)
2939 if (last_action == Insert)
2941 const int ierr =
matrix->GlobalAssemble(*column_space_map,
2953 add(i, 1, &j, &
value,
false);
2962# ifndef DEAL_II_WITH_64BIT_INDICES
2963 return matrix->NumGlobalRows();
2965 return matrix->NumGlobalRows64();
2984 return matrix->NumMyRows();
2989 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
2993# ifndef DEAL_II_WITH_64BIT_INDICES
2998 end =
matrix->RowMap().MaxMyGID64() + 1;
3012 return static_cast<std::uint64_t
>(
matrix->NumGlobalNonzeros64());
3017 template <
typename SparsityPatternType>
3019 const SparsityPatternType &sparsity_pattern,
3020 const MPI_Comm communicator,
3021 const bool exchange_data)
3023 reinit(parallel_partitioning,
3024 parallel_partitioning,
3032 template <
typename number>
3034 const IndexSet ¶llel_partitioning,
3035 const ::SparseMatrix<number> &sparse_matrix,
3036 const MPI_Comm communicator,
3037 const double drop_tolerance,
3038 const bool copy_values,
3039 const ::SparsityPattern *use_this_sparsity)
3043 reinit(parallel_partitioning,
3044 parallel_partitioning,
3055 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
3069 return IndexSet(
matrix->DomainMap());
3076 return IndexSet(
matrix->RangeMap());
3095 template <
typename VectorType>
3097 const VectorType &x,
3098 const VectorType &b)
const
3104 return dst.l2_norm();
3110 namespace LinearOperatorImplementation
3112 template <
typename EpetraOpType>
3113 TrilinosPayload::TrilinosPayload(
3115 const bool supports_inverse_operations,
3116 const bool use_transpose,
3117 const MPI_Comm mpi_communicator,
3118 const IndexSet &locally_owned_domain_indices,
3119 const IndexSet &locally_owned_range_indices)
3120 : use_transpose(use_transpose)
3121 , communicator(mpi_communicator)
3123 locally_owned_domain_indices.make_trilinos_map(communicator.Comm()))
3125 locally_owned_range_indices.make_trilinos_map(communicator.Comm()))
3130 Assert(&tril_src != &tril_dst,
3137 const int ierr = op.Apply(tril_src, tril_dst);
3144 Assert(&tril_src != &tril_dst,
3149 !op.UseTranspose());
3151 op.SetUseTranspose(!op.UseTranspose());
3152 const int ierr = op.Apply(tril_src, tril_dst);
3154 op.SetUseTranspose(!op.UseTranspose());
3157 if (supports_inverse_operations)
3163 &tril_src != &tril_dst,
3168 !op.UseTranspose());
3170 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3178 &tril_src != &tril_dst,
3185 op.SetUseTranspose(!op.UseTranspose());
3186 const int ierr = op.ApplyInverse(tril_src, tril_dst);
3188 op.SetUseTranspose(!op.UseTranspose());
3196 "Uninitialized TrilinosPayload::inv_vmult called. "
3197 "The operator does not support inverse operations."));
3203 "Uninitialized TrilinosPayload::inv_Tvmult called. "
3204 "The operator does not support inverse operations."));
3210 template <
typename Solver,
typename Preconditioner>
3212 std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
3213 std::is_base_of_v<TrilinosWrappers::PreconditionBase, Preconditioner>,
3217 const Preconditioner &preconditioner)
const
3219 const auto &payload = *
this;
3225 return_op.inv_vmult = [payload, &solver, &preconditioner](
3230 Assert(&tril_src != &tril_dst,
3235 !payload.UseTranspose());
3236 solver.solve(payload, tril_dst, tril_src, preconditioner);
3239 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3244 Assert(&tril_src != &tril_dst,
3249 payload.UseTranspose());
3252 solver.solve(payload, tril_dst, tril_src, preconditioner);
3258 if (return_op.UseTranspose() ==
true)
3259 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3264 template <
typename Solver,
typename Preconditioner>
3266 !(std::is_base_of_v<TrilinosWrappers::SolverBase, Solver> &&
3267 std::is_base_of_v<TrilinosWrappers::PreconditionBase,
3277 ExcMessage(
"Payload inv_vmult disabled because of "
3278 "incompatible solver/preconditioner choice."));
3284 ExcMessage(
"Payload inv_vmult disabled because of "
3285 "incompatible solver/preconditioner choice."));
3295 const size_type n_cols,
3296 const size_type *col_indices,
3298 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)
::types::global_dof_index size_type
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
const SparseMatrix MatrixType
TrilinosScalar value() const
const Accessor< Constness > & operator*() const
const Accessor< Constness > * operator->() const
::types::global_dof_index difference_type
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
TrilinosScalar residual(VectorType &dst, const VectorType &x, const VectorType &b) const
SparseMatrix & operator*=(const TrilinosScalar factor)
SparseMatrixIterators::Iterator< false > iterator
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
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > Tvmult(VectorType &dst, const VectorType &src) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
void clear_rows(const ArrayView< const size_type > &rows, 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
SparseMatrix & operator/=(const TrilinosScalar factor)
void Tvmult_add(VectorType &dst, const VectorType &src) const
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)
bool is_compressed() const
TrilinosScalar el(const size_type i, const size_type j) const
IndexSet locally_owned_domain_indices() 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
std::enable_if_t< std::is_same_v< typename VectorType::value_type, TrilinosScalar > > vmult(VectorType &dst, const VectorType &src) 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
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)
SparseMatrixIterators::Iterator< true > const_iterator
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
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
::types::global_dof_index size_type
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
std::enable_if_t< std::is_base_of_v< TrilinosWrappers::SolverBase, Solver > &&std::is_base_of_v< TrilinosWrappers::PreconditionBase, Preconditioner >, TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
TrilinosPayload identity_payload() 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
std::enable_if_t< !(std::is_base_of_v< TrilinosWrappers::SolverBase, Solver > &&std::is_base_of_v< TrilinosWrappers::PreconditionBase, Preconditioner >), TrilinosPayload > inverse_payload(Solver &, const Preconditioner &) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
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)
const bool IsBlockVector< VectorType >::value
types::global_dof_index size_type
@ matrix
Contents is actually a matrix.
SparseMatrix< double > SparseMatrix
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
VectorTraits::size_type 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 >::difference_type difference_type
typename ::TrilinosWrappers::SparseMatrixIterators::Iterator< Constness >::value_type value_type
forward_iterator_tag iterator_category