16 #ifndef dealii_trilinos_sparse_matrix_h 17 # define dealii_trilinos_sparse_matrix_h 20 # include <deal.II/base/config.h> 22 # ifdef DEAL_II_WITH_TRILINOS 24 # include <deal.II/base/index_set.h> 25 # include <deal.II/base/subscriptor.h> 27 # include <deal.II/lac/exceptions.h> 28 # include <deal.II/lac/full_matrix.h> 29 # include <deal.II/lac/trilinos_epetra_vector.h> 30 # include <deal.II/lac/trilinos_tpetra_vector.h> 31 # include <deal.II/lac/trilinos_vector.h> 32 # include <deal.II/lac/vector_memory.h> 33 # include <deal.II/lac/vector_operation.h> 35 # include <Epetra_Comm.h> 36 # include <Epetra_CrsGraph.h> 37 # include <Epetra_Export.h> 38 # include <Epetra_FECrsMatrix.h> 39 # include <Epetra_Map.h> 40 # include <Epetra_MultiVector.h> 41 # include <Epetra_Operator.h> 45 # include <type_traits> 47 # ifdef DEAL_II_WITH_MPI 48 # include <Epetra_MpiComm.h> 51 # include <Epetra_SerialComm.h> 54 DEAL_II_NAMESPACE_OPEN
57 template <
typename MatrixType>
60 template <
typename number>
79 template <
bool Constness>
94 <<
"You tried to access row " << arg1
95 <<
" of a distributed sparsity pattern, " 96 <<
" but only rows " << arg2 <<
" through " << arg3
97 <<
" are stored locally and can be accessed.");
198 template <
bool Constess>
237 template <
bool Other>
271 operator TrilinosScalar()
const;
277 operator=(
const TrilinosScalar n)
const;
283 operator+=(
const TrilinosScalar n)
const;
289 operator-=(
const TrilinosScalar n)
const;
295 operator*=(
const TrilinosScalar n)
const;
301 operator/=(
const TrilinosScalar n)
const;
339 friend class Reference;
356 template <
bool Constness>
380 template <
bool Other>
424 operator<(const Iterator<Constness> &)
const;
438 <<
"Attempt to access element " << arg2 <<
" of row " 439 << arg1 <<
" which doesn't have that many elements.");
447 template <
bool Other>
528 <<
"You tried to access row " << arg1
529 <<
" of a non-contiguous locally owned row set." 530 <<
" The row " << arg1
531 <<
" is not stored locally and can't be accessed.");
582 const unsigned int n_max_entries_per_row);
593 const std::vector<unsigned int> &n_entries_per_row);
637 template <
typename SparsityPatternType>
639 reinit(
const SparsityPatternType &sparsity_pattern);
687 template <
typename number>
689 reinit(const ::SparseMatrix<number> &dealii_sparse_matrix,
690 const double drop_tolerance = 1e-13,
691 const bool copy_values =
true,
692 const ::SparsityPattern * use_this_sparsity =
nullptr);
700 reinit(
const Epetra_CrsMatrix &input_matrix,
const bool copy_values =
true);
722 const size_type n_max_entries_per_row = 0);
735 const std::vector<unsigned int> &n_entries_per_row);
756 SparseMatrix(
const Epetra_Map &row_parallel_partitioning,
757 const Epetra_Map &col_parallel_partitioning,
758 const size_type n_max_entries_per_row = 0);
777 SparseMatrix(
const Epetra_Map & row_parallel_partitioning,
778 const Epetra_Map & col_parallel_partitioning,
779 const std::vector<unsigned int> &n_entries_per_row);
807 template <
typename SparsityPatternType>
808 DEAL_II_DEPRECATED
void 809 reinit(
const Epetra_Map & parallel_partitioning,
810 const SparsityPatternType &sparsity_pattern,
811 const bool exchange_data =
false);
827 template <
typename SparsityPatternType>
828 DEAL_II_DEPRECATED
void 829 reinit(
const Epetra_Map & row_parallel_partitioning,
830 const Epetra_Map & col_parallel_partitioning,
831 const SparsityPatternType &sparsity_pattern,
832 const bool exchange_data =
false);
852 template <
typename number>
853 DEAL_II_DEPRECATED
void 854 reinit(
const Epetra_Map & parallel_partitioning,
855 const ::SparseMatrix<number> &dealii_sparse_matrix,
856 const double drop_tolerance = 1e-13,
857 const bool copy_values =
true,
858 const ::SparsityPattern * use_this_sparsity =
nullptr);
875 template <
typename number>
876 DEAL_II_DEPRECATED
void 877 reinit(
const Epetra_Map & row_parallel_partitioning,
878 const Epetra_Map & col_parallel_partitioning,
879 const ::SparseMatrix<number> &dealii_sparse_matrix,
880 const double drop_tolerance = 1e-13,
881 const bool copy_values =
true,
882 const ::SparsityPattern * use_this_sparsity =
nullptr);
901 const MPI_Comm & communicator = MPI_COMM_WORLD,
902 const unsigned int n_max_entries_per_row = 0);
912 const MPI_Comm & communicator,
913 const std::vector<unsigned int> &n_entries_per_row);
930 const IndexSet &col_parallel_partitioning,
931 const MPI_Comm &communicator = MPI_COMM_WORLD,
932 const size_type n_max_entries_per_row = 0);
949 const IndexSet & col_parallel_partitioning,
950 const MPI_Comm & communicator,
951 const std::vector<unsigned int> &n_entries_per_row);
973 template <
typename SparsityPatternType>
976 const SparsityPatternType &sparsity_pattern,
977 const MPI_Comm & communicator = MPI_COMM_WORLD,
978 const bool exchange_data =
false);
992 template <
typename SparsityPatternType>
995 const IndexSet & col_parallel_partitioning,
996 const SparsityPatternType &sparsity_pattern,
997 const MPI_Comm & communicator = MPI_COMM_WORLD,
998 const bool exchange_data =
false);
1016 template <
typename number>
1019 const ::SparseMatrix<number> &dealii_sparse_matrix,
1020 const MPI_Comm & communicator = MPI_COMM_WORLD,
1021 const double drop_tolerance = 1e-13,
1022 const bool copy_values =
true,
1023 const ::SparsityPattern * use_this_sparsity =
nullptr);
1038 template <
typename number>
1041 const IndexSet & col_parallel_partitioning,
1042 const ::SparseMatrix<number> &dealii_sparse_matrix,
1043 const MPI_Comm & communicator = MPI_COMM_WORLD,
1044 const double drop_tolerance = 1e-13,
1045 const bool copy_values =
true,
1046 const ::SparsityPattern * use_this_sparsity =
nullptr);
1084 std::pair<size_type, size_type>
1245 set(
const std::vector<size_type> & indices,
1247 const bool elide_zero_values =
false);
1255 set(
const std::vector<size_type> & row_indices,
1256 const std::vector<size_type> & col_indices,
1258 const bool elide_zero_values =
false);
1289 const std::vector<size_type> & col_indices,
1290 const std::vector<TrilinosScalar> &values,
1291 const bool elide_zero_values =
false);
1320 template <
typename Number>
1325 const Number * values,
1326 const bool elide_zero_values =
false);
1359 add(
const std::vector<size_type> & indices,
1361 const bool elide_zero_values =
true);
1369 add(
const std::vector<size_type> & row_indices,
1370 const std::vector<size_type> & col_indices,
1372 const bool elide_zero_values =
true);
1389 const std::vector<size_type> & col_indices,
1390 const std::vector<TrilinosScalar> &values,
1391 const bool elide_zero_values =
true);
1410 const TrilinosScalar *values,
1411 const bool elide_zero_values =
true,
1412 const bool col_indices_are_sorted =
false);
1492 clear_rows(
const std::vector<size_type> &rows,
1493 const TrilinosScalar new_diag_value = 0);
1588 template <
typename VectorType>
1589 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1590 TrilinosScalar>::value>::type
1591 vmult(VectorType &dst,
const VectorType &src)
const;
1599 template <
typename VectorType>
1600 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1601 TrilinosScalar>::value>::type
1602 vmult(VectorType &dst,
const VectorType &src)
const;
1618 template <
typename VectorType>
1619 typename std::enable_if<std::is_same<
typename VectorType::value_type,
1620 TrilinosScalar>::value>::type
1621 Tvmult(VectorType &dst,
const VectorType &src)
const;
1629 template <
typename VectorType>
1630 typename std::enable_if<!std::is_same<
typename VectorType::value_type,
1631 TrilinosScalar>::value>::type
1632 Tvmult(VectorType &dst,
const VectorType &src)
const;
1643 template <
typename VectorType>
1645 vmult_add(VectorType &dst,
const VectorType &src)
const;
1657 template <
typename VectorType>
1659 Tvmult_add(VectorType &dst,
const VectorType &src)
const;
1813 const Epetra_CrsMatrix &
1820 const Epetra_CrsGraph &
2015 print(std::ostream &out,
2016 const bool write_extended_trilinos_info =
false)
const;
2029 <<
"An error with error number " << arg1
2030 <<
" occurred while calling a Trilinos function");
2038 <<
"The entry with index <" << arg1 <<
',' << arg2
2039 <<
"> does not exist.");
2045 "You are attempting an operation on two matrices that " 2046 "are the same object, but the operation requires that the " 2047 "two objects are in fact different.");
2062 <<
"You tried to access element (" << arg1 <<
"/" << arg2
2064 <<
" of a distributed matrix, but only rows " << arg3
2065 <<
" through " << arg4
2066 <<
" are stored locally and can be accessed.");
2074 <<
"You tried to access element (" << arg1 <<
"/" << arg2
2076 <<
" of a sparse matrix, but it appears to not" 2077 <<
" exist in the Trilinos sparsity pattern.");
2166 check_vector_map_equality(
const Epetra_CrsMatrix & mtrx,
2167 const Epetra_MultiVector &src,
2168 const Epetra_MultiVector &dst,
2173 Assert(src.Map().SameAs(mtrx.DomainMap()) ==
true,
2175 "Column map of matrix does not fit with vector map!"));
2176 Assert(dst.Map().SameAs(mtrx.RangeMap()) ==
true,
2177 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2181 Assert(src.Map().SameAs(mtrx.RangeMap()) ==
true,
2183 "Column map of matrix does not fit with vector map!"));
2184 Assert(dst.Map().SameAs(mtrx.DomainMap()) ==
true,
2185 ExcMessage(
"Row map of matrix does not fit with vector map!"));
2193 check_vector_map_equality(
const Epetra_Operator & op,
2194 const Epetra_MultiVector &src,
2195 const Epetra_MultiVector &dst,
2200 Assert(src.Map().SameAs(op.OperatorDomainMap()) ==
true,
2202 "Column map of operator does not fit with vector map!"));
2203 Assert(dst.Map().SameAs(op.OperatorRangeMap()) ==
true,
2205 "Row map of operator does not fit with vector map!"));
2209 Assert(src.Map().SameAs(op.OperatorRangeMap()) ==
true,
2211 "Column map of operator does not fit with vector map!"));
2212 Assert(dst.Map().SameAs(op.OperatorDomainMap()) ==
true,
2214 "Row map of operator does not fit with vector map!"));
2221 namespace LinearOperatorImplementation
2349 template <
typename Solver,
typename Preconditioner>
2350 typename std::enable_if<
2351 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2353 Preconditioner>::value,
2374 template <
typename Solver,
typename Preconditioner>
2375 typename std::enable_if<
2376 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
2378 Preconditioner>::value),
2427 std::function<void(VectorType &, const VectorType &)>
vmult;
2436 std::function<void(VectorType &, const VectorType &)>
Tvmult;
2446 std::function<void(VectorType &, const VectorType &)>
inv_vmult;
2456 std::function<void(VectorType &, const VectorType &)>
inv_Tvmult;
2539 virtual const char *
2540 Label()
const override;
2549 virtual const Epetra_Comm &
2550 Comm()
const override;
2559 virtual const Epetra_Map &
2570 virtual const Epetra_Map &
2585 # ifdef DEAL_II_WITH_MPI 2658 visit_present_row();
2662 inline AccessorBase::size_type
2663 AccessorBase::row()
const 2670 inline AccessorBase::size_type
2671 AccessorBase::column()
const 2674 return (*colnum_cache)[a_index];
2678 inline AccessorBase::size_type
2679 AccessorBase::index()
const 2686 inline Accessor<true>::Accessor(MatrixType * matrix,
2687 const size_type row,
2688 const size_type index)
2693 template <
bool Other>
2694 inline Accessor<true>::Accessor(
const Accessor<Other> &other)
2695 : AccessorBase(other)
2699 inline TrilinosScalar
2700 Accessor<true>::value()
const 2703 return (*value_cache)[a_index];
2707 inline Accessor<false>::Reference::Reference(
const Accessor<false> &acc)
2708 : accessor(const_cast<Accessor<false> &>(acc))
2712 inline Accessor<false>::Reference::operator TrilinosScalar()
const 2714 return (*accessor.value_cache)[accessor.a_index];
2717 inline const Accessor<false>::Reference &
2718 Accessor<false>::Reference::operator=(
const TrilinosScalar n)
const 2720 (*accessor.value_cache)[accessor.a_index] = n;
2721 accessor.matrix->set(accessor.row(),
2723 static_cast<TrilinosScalar
>(*this));
2728 inline const Accessor<false>::Reference &
2729 Accessor<false>::Reference::operator+=(
const TrilinosScalar n)
const 2731 (*accessor.value_cache)[accessor.a_index] += n;
2732 accessor.matrix->set(accessor.row(),
2734 static_cast<TrilinosScalar
>(*this));
2739 inline const Accessor<false>::Reference &
2740 Accessor<false>::Reference::operator-=(
const TrilinosScalar n)
const 2742 (*accessor.value_cache)[accessor.a_index] -= n;
2743 accessor.matrix->set(accessor.row(),
2745 static_cast<TrilinosScalar
>(*this));
2750 inline const Accessor<false>::Reference &
2751 Accessor<false>::Reference::operator*=(
const TrilinosScalar n)
const 2753 (*accessor.value_cache)[accessor.a_index] *= n;
2754 accessor.matrix->set(accessor.row(),
2756 static_cast<TrilinosScalar
>(*this));
2761 inline const Accessor<false>::Reference &
2762 Accessor<false>::Reference::operator/=(
const TrilinosScalar n)
const 2764 (*accessor.value_cache)[accessor.a_index] /= n;
2765 accessor.matrix->set(accessor.row(),
2767 static_cast<TrilinosScalar
>(*this));
2772 inline Accessor<false>::Accessor(MatrixType * matrix,
2773 const size_type row,
2774 const size_type index)
2775 : AccessorBase(
matrix, row, index)
2779 inline Accessor<false>::Reference
2780 Accessor<false>::value()
const 2788 template <
bool Constness>
2789 inline Iterator<Constness>::Iterator(MatrixType * matrix,
2790 const size_type row,
2791 const size_type index)
2792 : accessor(
matrix, row, index)
2796 template <
bool Constness>
2797 template <
bool Other>
2798 inline Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2799 : accessor(other.accessor)
2803 template <
bool Constness>
2804 inline Iterator<Constness> &
2805 Iterator<Constness>::operator++()
2815 if (accessor.a_index >= accessor.colnum_cache->size())
2817 accessor.a_index = 0;
2820 while ((accessor.a_row < accessor.matrix->m()) &&
2821 ((accessor.matrix->in_local_range(accessor.a_row) ==
false) ||
2822 (accessor.matrix->row_length(accessor.a_row) == 0)))
2825 accessor.visit_present_row();
2831 template <
bool Constness>
2832 inline Iterator<Constness>
2833 Iterator<Constness>::operator++(
int)
2835 const Iterator<Constness> old_state = *
this;
2842 template <
bool Constness>
2850 template <
bool Constness>
2851 inline const Accessor<Constness> *Iterator<Constness>::operator->()
const 2858 template <
bool Constness>
2860 Iterator<Constness>::operator==(
const Iterator<Constness> &other)
const 2862 return (accessor.a_row == other.accessor.a_row &&
2863 accessor.a_index == other.accessor.a_index);
2868 template <
bool Constness>
2870 Iterator<Constness>::operator!=(
const Iterator<Constness> &other)
const 2872 return !(*
this == other);
2877 template <
bool Constness>
2879 Iterator<Constness>::operator<(
const Iterator<Constness> &other)
const 2881 return (accessor.row() < other.accessor.row() ||
2882 (accessor.row() == other.accessor.row() &&
2883 accessor.index() < other.accessor.index()));
2887 template <
bool Constness>
2889 Iterator<Constness>::operator>(
const Iterator<Constness> &other)
const 2891 return (other < *
this);
2909 return const_iterator(
this, m(), 0);
2918 if (in_local_range(r) && (row_length(r) > 0))
2919 return const_iterator(
this, r, 0);
2934 for (size_type i = r + 1; i < m(); ++i)
2935 if (in_local_range(i) && (row_length(i) > 0))
2936 return const_iterator(
this, i, 0);
2956 return iterator(
this, m(), 0);
2965 if (in_local_range(r) && (row_length(r) > 0))
2966 return iterator(
this, r, 0);
2981 for (size_type i = r + 1; i < m(); ++i)
2982 if (in_local_range(i) && (row_length(i) > 0))
2983 return iterator(
this, i, 0);
2995 TrilinosWrappers::types::int_type begin, end;
2996 # ifndef DEAL_II_WITH_64BIT_INDICES 2997 begin =
matrix->RowMap().MinMyGID();
2998 end =
matrix->RowMap().MaxMyGID() + 1;
3000 begin =
matrix->RowMap().MinMyGID64();
3001 end =
matrix->RowMap().MaxMyGID64() + 1;
3004 return ((index >= static_cast<size_type>(begin)) &&
3005 (index < static_cast<size_type>(end)));
3023 SparseMatrix::set<TrilinosScalar>(
const size_type row,
3024 const size_type n_cols,
3025 const size_type * col_indices,
3026 const TrilinosScalar *values,
3027 const bool elide_zero_values);
3031 template <
typename Number>
3034 const size_type n_cols,
3035 const size_type *col_indices,
3036 const Number * values,
3037 const bool elide_zero_values)
3039 std::vector<TrilinosScalar> trilinos_values(n_cols);
3040 std::copy(values, values + n_cols, trilinos_values.begin());
3042 row, n_cols, col_indices, trilinos_values.data(), elide_zero_values);
3050 const TrilinosScalar value)
3054 set(i, 1, &j, &value,
false);
3062 const bool elide_zero_values)
3064 Assert(indices.size() == values.
m(),
3068 for (size_type i = 0; i < indices.size(); ++i)
3081 const TrilinosScalar value)
3093 if (last_action == Insert)
3096 ierr =
matrix->GlobalAssemble(*column_space_map,
3109 add(i, 1, &j, &value,
false);
3119 # ifndef DEAL_II_WITH_64BIT_INDICES 3120 return matrix->NumGlobalRows();
3122 return matrix->NumGlobalRows64();
3135 # ifndef DEAL_II_WITH_64BIT_INDICES 3136 return column_space_map->NumGlobalElements();
3138 return column_space_map->NumGlobalElements64();
3147 return matrix->NumMyRows();
3152 inline std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
3155 size_type begin, end;
3156 # ifndef DEAL_II_WITH_64BIT_INDICES 3157 begin =
matrix->RowMap().MinMyGID();
3158 end =
matrix->RowMap().MaxMyGID() + 1;
3160 begin =
matrix->RowMap().MinMyGID64();
3161 end =
matrix->RowMap().MaxMyGID64() + 1;
3164 return std::make_pair(begin, end);
3172 # ifndef DEAL_II_WITH_64BIT_INDICES 3173 return matrix->NumGlobalNonzeros();
3175 return matrix->NumGlobalNonzeros64();
3181 template <
typename SparsityPatternType>
3184 const SparsityPatternType &sparsity_pattern,
3185 const MPI_Comm & communicator,
3186 const bool exchange_data)
3188 reinit(parallel_partitioning,
3189 parallel_partitioning,
3197 template <
typename number>
3200 const ::SparseMatrix<number> &sparse_matrix,
3201 const MPI_Comm & communicator,
3202 const double drop_tolerance,
3203 const bool copy_values,
3204 const ::SparsityPattern * use_this_sparsity)
3208 reinit(parallel_partitioning,
3209 parallel_partitioning,
3218 inline const Epetra_CrsMatrix &
3221 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
3226 inline const Epetra_CrsGraph &
3267 namespace LinearOperatorImplementation
3269 template <
typename Solver,
typename Preconditioner>
3270 typename std::enable_if<
3271 std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3273 Preconditioner>::value,
3274 TrilinosPayload>::type
3275 TrilinosPayload::inverse_payload(
3277 const Preconditioner &preconditioner)
const 3279 const auto &payload = *
this;
3281 TrilinosPayload return_op(payload);
3285 return_op.inv_vmult = [payload, &solver, &preconditioner](
3286 TrilinosPayload::Domain & tril_dst,
3287 const TrilinosPayload::Range &tril_src) {
3290 Assert(&tril_src != &tril_dst,
3292 internal::check_vector_map_equality(payload,
3295 !payload.UseTranspose());
3296 solver.solve(payload, tril_dst, tril_src, preconditioner);
3299 return_op.inv_Tvmult = [payload, &solver, &preconditioner](
3300 TrilinosPayload::Range & tril_dst,
3301 const TrilinosPayload::Domain &tril_src) {
3304 Assert(&tril_src != &tril_dst,
3306 internal::check_vector_map_equality(payload,
3309 payload.UseTranspose());
3311 const_cast<TrilinosPayload &
>(payload).
transpose();
3312 solver.solve(payload, tril_dst, tril_src, preconditioner);
3313 const_cast<TrilinosPayload &
>(payload).
transpose();
3318 if (return_op.UseTranspose() ==
true)
3319 std::swap(return_op.inv_vmult, return_op.inv_Tvmult);
3324 template <
typename Solver,
typename Preconditioner>
3325 typename std::enable_if<
3326 !(std::is_base_of<TrilinosWrappers::SolverBase, Solver>::value &&
3328 Preconditioner>::value),
3329 TrilinosPayload>::type
3330 TrilinosPayload::inverse_payload(
Solver &,
const Preconditioner &)
const 3332 TrilinosPayload return_op(*
this);
3334 return_op.inv_vmult = [](TrilinosPayload::Domain &,
3335 const TrilinosPayload::Range &) {
3337 ExcMessage(
"Payload inv_vmult disabled because of " 3338 "incompatible solver/preconditioner choice."));
3341 return_op.inv_Tvmult = [](TrilinosPayload::Range &,
3342 const TrilinosPayload::Domain &) {
3344 ExcMessage(
"Payload inv_vmult disabled because of " 3345 "incompatible solver/preconditioner choice."));
3355 SparseMatrix::set<TrilinosScalar>(
const size_type row,
3356 const size_type n_cols,
3357 const size_type * col_indices,
3358 const TrilinosScalar *values,
3359 const bool elide_zero_values);
3365 DEAL_II_NAMESPACE_CLOSE
3368 # endif // DEAL_II_WITH_TRILINOS static const bool zero_addition_can_be_elided
TrilinosScalar diag_element(const size_type i) const
static ::ExceptionBase & ExcTrilinosError(int arg1)
TrilinosScalar el(const size_type i, const size_type j) const
virtual double NormInf() const override
const_iterator end() const
const Epetra_Map & col_partitioner() const
bool operator>(const Iterator< Constness > &) const
TrilinosScalar frobenius_norm() const
void mmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
#define DeclException2(Exception2, type1, type2, outsequence)
TrilinosPayload null_payload() const
virtual int SetUseTranspose(bool UseTranspose) override
Contents is actually a matrix.
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcBeyondEndOfMatrix()
bool operator==(const Iterator< Constness > &) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
const Epetra_Map & range_partitioner() const
void vmult_add(VectorType &dst, const VectorType &src) const
Epetra_MultiVector VectorType
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
Epetra_CombineMode last_action
bool is_compressed() const
virtual const char * Label() const override
TrilinosScalar operator()(const size_type i, const size_type j) const
std::function< void(VectorType &, const VectorType &)> Tvmult
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1, std::size_t arg2, std::size_t arg3)
const Epetra_Map & row_partitioner() const
#define AssertThrow(cond, exc)
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
SymmetricTensor< rank_, dim, typename ProductType< Number, OtherNumber >::type > operator+(const SymmetricTensor< rank_, dim, Number > &left, const SymmetricTensor< rank_, dim, OtherNumber > &right)
TrilinosScalar residual(MPI::Vector &dst, const MPI::Vector &x, const MPI::Vector &b) const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcTrilinosError(int arg1)
SymmetricTensor< rank_, dim, Number > operator*(const SymmetricTensor< rank_, dim, Number > &t, const Number &factor)
TrilinosScalar matrix_scalar_product(const MPI::Vector &u, const MPI::Vector &v) const
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
std::shared_ptr< std::vector< size_type > > colnum_cache
SparseMatrix & operator*=(const TrilinosScalar factor)
std::unique_ptr< Epetra_FECrsMatrix > matrix
virtual int Apply(const VectorType &X, VectorType &Y) const override
SparseMatrix & operator=(const SparseMatrix &)=delete
::types::global_dof_index size_type
TrilinosPayload transpose_payload() const
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1)
::types::global_dof_index size_type
bool in_local_range(const size_type index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
virtual const Epetra_Map & OperatorRangeMap() const override
const Accessor< Constness > & operator*() const
TrilinosPayload identity_payload() const
#define DeclException1(Exception1, type1, outsequence)
void add(const size_type i, const size_type j, const TrilinosScalar value)
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
std::unique_ptr< Epetra_Export > nonlocal_matrix_exporter
size_type memory_consumption() const
#define DeclExceptionMsg(Exception, defaulttext)
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
bool operator!=(const Iterator< Constness > &) const
Epetra_MpiComm communicator
#define DeclException0(Exception0)
typename Accessor< Constness >::MatrixType MatrixType
std::function< void(VectorType &, const VectorType &)> inv_vmult
TrilinosScalar matrix_norm_square(const MPI::Vector &v) const
const Accessor< Constness > * operator->() const
virtual int ApplyInverse(const VectorType &Y, VectorType &X) 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_range_indices() const
const Epetra_CrsMatrix & trilinos_matrix() const
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
virtual bool UseTranspose() const override
AccessorBase(SparseMatrix *matrix, const size_type row, const size_type index)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
virtual ~SparseMatrix() override=default
IndexSet locally_owned_range_indices() const
virtual const Epetra_Map & OperatorDomainMap() const override
size_type n_nonzero_elements() const
std::enable_if< std::is_same< typename VectorType::value_type, TrilinosScalar >::value >::type Tvmult(VectorType &dst, const VectorType &src) const
MPI_Comm get_mpi_communicator() const
static ::ExceptionBase & ExcIteratorPastEnd()
types::global_dof_index size_type
const_iterator begin() const
static ::ExceptionBase & ExcNotQuadratic()
void swap(Vector< Number > &u, Vector< Number > &v)
unsigned int row_length(const size_type row) const
Accessor< Constness > accessor
void Tvmult_add(VectorType &dst, const VectorType &src) const
unsigned int global_dof_index
virtual ~TrilinosPayload() override=default
std::enable_if< std::is_same< typename VectorType::value_type, TrilinosScalar >::value >::type vmult(VectorType &dst, const VectorType &src) const
TrilinosScalar value_type
void reinit(const SparsityPatternType &sparsity_pattern)
TrilinosScalar linfty_norm() const
Iterator< Constness > & operator++()
void copy_from(const SparseMatrix &source)
const Epetra_Map & domain_partitioner() const
std::function< void(VectorType &, const VectorType &)> vmult
SparseMatrix & operator/=(const TrilinosScalar factor)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
virtual const Epetra_Comm & Comm() const override
std::pair< size_type, size_type > local_range() const
IndexSet locally_owned_domain_indices() const
static ::ExceptionBase & ExcMatrixNotCompressed()
IndexSet locally_owned_domain_indices() const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const MPI::Vector &V=MPI::Vector()) const
#define DeclException3(Exception3, type1, type2, type3, outsequence)
TrilinosScalar value() const
std::unique_ptr< Epetra_Map > column_space_map
virtual bool HasNormInf() const override
void compress(::VectorOperation::values operation)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
std::unique_ptr< Epetra_CrsMatrix > nonlocal_matrix
MPI_Comm get_mpi_communicator() const
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
unsigned int local_size() const
TrilinosScalar l1_norm() const
::types::global_dof_index size_type
static ::ExceptionBase & ExcInternalError()
std::shared_ptr< std::vector< TrilinosScalar > > value_cache