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/std_cxx11/shared_ptr.h> 25 # include <deal.II/base/subscriptor.h> 26 # include <deal.II/base/index_set.h> 27 # include <deal.II/lac/full_matrix.h> 28 # include <deal.II/lac/exceptions.h> 29 # include <deal.II/lac/trilinos_vector.h> 30 # include <deal.II/lac/vector_view.h> 36 # define TrilinosScalar double 38 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
39 # include <Epetra_FECrsMatrix.h> 40 # include <Epetra_Map.h> 41 # include <Epetra_CrsGraph.h> 42 # include <Epetra_MultiVector.h> 43 # ifdef DEAL_II_WITH_MPI 44 # include <Epetra_MpiComm.h> 47 # include "Epetra_SerialComm.h" 53 DEAL_II_NAMESPACE_OPEN
86 std::size_t, std::size_t, std::size_t,
87 <<
"You tried to access row " << arg1
88 <<
" of a distributed sparsity pattern, " 89 <<
" but only rows " << arg2 <<
" through " << arg3
90 <<
" are stored locally and can be accessed.");
115 const size_type
index);
120 size_type
row()
const;
125 size_type
index()
const;
187 template <
bool Constess>
193 TrilinosScalar value()
const;
198 TrilinosScalar &value();
226 template <
bool Other>
232 TrilinosScalar value()
const;
258 operator TrilinosScalar ()
const;
263 const Reference &operator = (
const TrilinosScalar n)
const;
268 const Reference &operator += (
const TrilinosScalar n)
const;
273 const Reference &operator -= (
const TrilinosScalar n)
const;
278 const Reference &operator *= (
const TrilinosScalar n)
const;
283 const Reference &operator /= (
const TrilinosScalar n)
const;
311 Reference value()
const;
321 friend class Reference;
338 template <
bool Constness>
359 const size_type
index);
364 template <
bool Other>
403 bool operator < (const Iterator<Constness> &)
const;
414 size_type, size_type,
415 <<
"Attempt to access element " << arg2
416 <<
" of row " << arg1
417 <<
" which doesn't have that many elements.");
425 template <
bool Other>
friend class Iterator;
513 static const bool zero_addition_can_be_elided =
true;
549 const unsigned int n_max_entries_per_row);
560 const std::vector<unsigned int> &n_entries_per_row);
587 template<
typename SparsityPatternType>
588 void reinit (
const SparsityPatternType &sparsity_pattern);
634 template <
typename number>
635 void reinit (const ::SparseMatrix<number> &dealii_sparse_matrix,
636 const double drop_tolerance=1e-13,
637 const bool copy_values=
true,
638 const ::SparsityPattern *use_this_sparsity=0);
645 void reinit (
const Epetra_CrsMatrix &input_matrix,
646 const bool copy_values =
true);
699 SparseMatrix (const Epetra_Map &row_parallel_partitioning,
700 const Epetra_Map &col_parallel_partitioning,
701 const size_type n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
719 SparseMatrix (const Epetra_Map &row_parallel_partitioning,
720 const Epetra_Map &col_parallel_partitioning,
721 const
std::vector<
unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
749 template<typename SparsityPatternType>
750 void reinit (const Epetra_Map ¶llel_partitioning,
751 const SparsityPatternType &sparsity_pattern,
752 const
bool exchange_data = false) DEAL_II_DEPRECATED;
768 template<typename SparsityPatternType>
769 void reinit (const Epetra_Map &row_parallel_partitioning,
770 const Epetra_Map &col_parallel_partitioning,
771 const SparsityPatternType &sparsity_pattern,
772 const
bool exchange_data = false) DEAL_II_DEPRECATED;
792 template <typename number>
793 void reinit (const Epetra_Map ¶llel_partitioning,
795 const
double drop_tolerance=1e-13,
796 const
bool copy_values=true,
814 template <typename number>
815 void reinit (const Epetra_Map &row_parallel_partitioning,
816 const Epetra_Map &col_parallel_partitioning,
818 const
double drop_tolerance=1e-13,
819 const
bool copy_values=true,
839 const MPI_Comm &communicator = MPI_COMM_WORLD,
840 const
unsigned int n_max_entries_per_row = 0);
850 const MPI_Comm &communicator,
851 const
std::vector<
unsigned int> &n_entries_per_row);
868 const
IndexSet &col_parallel_partitioning,
869 const MPI_Comm &communicator = MPI_COMM_WORLD,
870 const size_type n_max_entries_per_row = 0);
887 const
IndexSet &col_parallel_partitioning,
888 const MPI_Comm &communicator,
889 const
std::vector<
unsigned int> &n_entries_per_row);
911 template<typename SparsityPatternType>
912 void reinit (const
IndexSet ¶llel_partitioning,
913 const SparsityPatternType &sparsity_pattern,
914 const MPI_Comm &communicator = MPI_COMM_WORLD,
915 const
bool exchange_data = false);
929 template<typename SparsityPatternType>
930 void reinit (const
IndexSet &row_parallel_partitioning,
931 const
IndexSet &col_parallel_partitioning,
932 const SparsityPatternType &sparsity_pattern,
933 const MPI_Comm &communicator = MPI_COMM_WORLD,
934 const
bool exchange_data = false);
952 template <typename number>
953 void reinit (const
IndexSet ¶llel_partitioning,
955 const MPI_Comm &communicator = MPI_COMM_WORLD,
956 const
double drop_tolerance=1e-13,
957 const
bool copy_values=true,
973 template <typename number>
974 void reinit (const
IndexSet &row_parallel_partitioning,
975 const
IndexSet &col_parallel_partitioning,
977 const MPI_Comm &communicator = MPI_COMM_WORLD,
978 const
double drop_tolerance=1e-13,
979 const
bool copy_values=true,
990 size_type m () const;
995 size_type n () const;
1005 unsigned int local_size () const;
1015 std::pair<size_type, size_type>
1016 local_range () const;
1022 bool in_local_range (const size_type
index) const;
1027 size_type n_nonzero_elements () const;
1032 unsigned int row_length (const size_type
row) const;
1040 bool is_compressed () const;
1047 size_type memory_consumption () const;
1052 MPI_Comm get_mpi_communicator () const;
1070 operator = (const
double d);
1131 void set (const size_type i,
1133 const TrilinosScalar value);
1167 void set (const
std::vector<size_type> &indices,
1168 const
FullMatrix<TrilinosScalar> &full_matrix,
1169 const
bool elide_zero_values = false);
1176 void set (const
std::vector<size_type> &row_indices,
1177 const
std::vector<size_type> &col_indices,
1178 const
FullMatrix<TrilinosScalar> &full_matrix,
1179 const
bool elide_zero_values = false);
1208 void set (const size_type
row,
1209 const
std::vector<size_type> &col_indices,
1210 const
std::vector<TrilinosScalar> &values,
1211 const
bool elide_zero_values = false);
1240 void set (const size_type row,
1241 const size_type n_cols,
1242 const size_type *col_indices,
1243 const TrilinosScalar *values,
1244 const
bool elide_zero_values = false);
1255 void add (const size_type i,
1257 const TrilinosScalar value);
1277 void add (const
std::vector<size_type> &indices,
1278 const
FullMatrix<TrilinosScalar> &full_matrix,
1279 const
bool elide_zero_values = true);
1286 void add (const
std::vector<size_type> &row_indices,
1287 const
std::vector<size_type> &col_indices,
1288 const
FullMatrix<TrilinosScalar> &full_matrix,
1289 const
bool elide_zero_values = true);
1304 void add (const size_type row,
1305 const
std::vector<size_type> &col_indices,
1306 const
std::vector<TrilinosScalar> &values,
1307 const
bool elide_zero_values = true);
1322 void add (const size_type row,
1323 const size_type n_cols,
1324 const size_type *col_indices,
1325 const TrilinosScalar *values,
1326 const
bool elide_zero_values = true,
1327 const
bool col_indices_are_sorted = false);
1332 SparseMatrix &operator *= (const TrilinosScalar factor);
1337 SparseMatrix &operator /= (const TrilinosScalar factor);
1351 void add (const TrilinosScalar factor,
1380 void clear_row (const size_type row,
1381 const TrilinosScalar new_diag_value = 0);
1403 void clear_rows (const
std::vector<size_type> &rows,
1404 const TrilinosScalar new_diag_value = 0);
1428 TrilinosScalar operator () (const size_type i,
1429 const size_type j) const;
1447 TrilinosScalar el (const size_type i,
1448 const size_type j) const;
1456 TrilinosScalar diag_element (const size_type i) const;
1485 template<typename VectorType>
1486 void vmult (VectorType &dst,
1487 const VectorType &src) const;
1511 template <typename VectorType>
1512 void Tvmult (VectorType &dst,
1513 const VectorType &src) const;
1538 template<typename VectorType>
1539 void vmult_add (VectorType &dst,
1540 const VectorType &src) const;
1565 template <typename VectorType>
1566 void Tvmult_add (VectorType &dst,
1567 const VectorType &src) const;
1593 TrilinosScalar matrix_norm_square (const
VectorBase &v) const;
1613 TrilinosScalar matrix_scalar_product (const
VectorBase &u,
1689 TrilinosScalar l1_norm () const;
1699 TrilinosScalar linfty_norm () const;
1705 TrilinosScalar frobenius_norm () const;
1717 const Epetra_CrsMatrix &trilinos_matrix () const;
1723 const Epetra_CrsGraph &trilinos_sparsity_pattern () const;
1732 const Epetra_Map &domain_partitioner () const DEAL_II_DEPRECATED;
1742 const Epetra_Map &range_partitioner () const DEAL_II_DEPRECATED;
1751 const Epetra_Map &row_partitioner () const DEAL_II_DEPRECATED;
1762 const Epetra_Map &col_partitioner () const DEAL_II_DEPRECATED;
1774 IndexSet locally_owned_domain_indices() const;
1781 IndexSet locally_owned_range_indices() const;
1808 const_iterator begin () const;
1819 const_iterator end () const;
1854 const_iterator begin (const size_type r) const;
1859 iterator begin (const size_type r);
1870 const_iterator end (const size_type r) const;
1875 iterator end (const size_type r);
1888 void write_ascii ();
1897 void print (
std::ostream &out,
1898 const
bool write_extended_trilinos_info = false) const;
1911 << "An error with error number " << arg1
1912 << " occurred while calling a Trilinos function");
1918 size_type, size_type,
1919 << "The entry with
index <" << arg1 << ',' << arg2
1920 << "> does not exist.");
1936 size_type, size_type, size_type, size_type,
1937 << "You tried to access element (" << arg1
1938 << "/" << arg2 << ")"
1939 << " of a distributed matrix, but only rows "
1940 << arg3 << " through " << arg4
1941 << " are stored locally and can be accessed.");
1947 size_type, size_type,
1948 << "You tried to access element (" << arg1
1949 << "/" << arg2 << ")"
1950 << " of a sparse matrix, but it appears to not"
1951 << " exist in the Trilinos sparsity pattern.");
1987 SparseMatrix &operator = (const SparseMatrix &);
2012 std_cxx11::shared_ptr<Epetra_Export> nonlocal_matrix_exporter;
2025 Epetra_CombineMode last_action;
2063 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2072 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2081 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2088 const size_type row,
2089 const size_type index)
2091 AccessorBase(const_cast<SparseMatrix *>(matrix), row, index)
2095 template <
bool Other>
2107 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2113 Accessor<false>::Reference::Reference (
2114 const Accessor<false> &acc)
2116 accessor(
const_cast<Accessor<false>&
>(acc))
2121 Accessor<false>::Reference::operator TrilinosScalar ()
const 2123 return (*accessor.value_cache)[accessor.a_index];
2128 Accessor<false>::Reference::operator = (
const TrilinosScalar n)
const 2130 (*accessor.value_cache)[accessor.a_index] = n;
2131 accessor.matrix->set(accessor.row(), accessor.column(),
2132 static_cast<TrilinosScalar
>(*this));
2139 Accessor<false>::Reference::operator += (
const TrilinosScalar n)
const 2141 (*accessor.value_cache)[accessor.a_index] += n;
2142 accessor.matrix->set(accessor.row(), accessor.column(),
2143 static_cast<TrilinosScalar
>(*this));
2150 Accessor<false>::Reference::operator -= (
const TrilinosScalar n)
const 2152 (*accessor.value_cache)[accessor.a_index] -= n;
2153 accessor.matrix->set(accessor.row(), accessor.column(),
2154 static_cast<TrilinosScalar
>(*this));
2161 Accessor<false>::Reference::operator *= (
const TrilinosScalar n)
const 2163 (*accessor.value_cache)[accessor.a_index] *= n;
2164 accessor.matrix->set(accessor.row(), accessor.column(),
2165 static_cast<TrilinosScalar
>(*this));
2172 Accessor<false>::Reference::operator /= (
const TrilinosScalar n)
const 2174 (*accessor.value_cache)[accessor.a_index] /= n;
2175 accessor.matrix->set(accessor.row(), accessor.column(),
2176 static_cast<TrilinosScalar
>(*this));
2183 const size_type row,
2184 const size_type index)
2194 Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2195 return Reference(*
this);
2200 template <
bool Constness>
2202 Iterator<Constness>::Iterator(MatrixType *matrix,
2203 const size_type row,
2204 const size_type index)
2206 accessor(matrix, row, index)
2210 template <
bool Constness>
2211 template <
bool Other>
2213 Iterator<Constness>::Iterator(
const Iterator<Other> &other)
2215 accessor(other.accessor)
2219 template <
bool Constness>
2221 Iterator<Constness> &
2224 Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
2232 if (accessor.a_index >= accessor.colnum_cache->size())
2234 accessor.a_index = 0;
2237 while ((accessor.a_row < accessor.matrix->m())
2239 ((accessor.matrix->in_local_range (accessor.a_row) ==
false)
2241 (accessor.matrix->row_length(accessor.a_row) == 0)))
2244 accessor.visit_present_row();
2250 template <
bool Constness>
2255 const Iterator<Constness> old_state = *
this;
2262 template <
bool Constness>
2264 const Accessor<Constness> &
2272 template <
bool Constness>
2274 const Accessor<Constness> *
2282 template <
bool Constness>
2287 return (accessor.a_row == other.accessor.a_row &&
2288 accessor.a_index == other.accessor.a_index);
2293 template <
bool Constness>
2298 return ! (*
this == other);
2303 template <
bool Constness>
2308 return (accessor.row() < other.accessor.row() ||
2309 (accessor.row() == other.accessor.row() &&
2310 accessor.index() < other.accessor.index()));
2314 template <
bool Constness>
2319 return (other < *
this);
2339 return const_iterator(
this, m(), 0);
2348 Assert (r < m(), ExcIndexRange(r, 0, m()));
2349 if (in_local_range (r)
2351 (row_length(r) > 0))
2352 return const_iterator(
this, r, 0);
2363 Assert (r < m(), ExcIndexRange(r, 0, m()));
2368 for (size_type i=r+1; i<m(); ++i)
2369 if (in_local_range (i)
2371 (row_length(i) > 0))
2372 return const_iterator(
this, i, 0);
2394 return iterator(
this, m(), 0);
2403 Assert (r < m(), ExcIndexRange(r, 0, m()));
2404 if (in_local_range (r)
2406 (row_length(r) > 0))
2407 return iterator(
this, r, 0);
2418 Assert (r < m(), ExcIndexRange(r, 0, m()));
2423 for (size_type i=r+1; i<m(); ++i)
2424 if (in_local_range (i)
2426 (row_length(i) > 0))
2427 return iterator(
this, i, 0);
2440 TrilinosWrappers::types::int_type begin, end;
2441 #ifndef DEAL_II_WITH_64BIT_INDICES 2442 begin = matrix->RowMap().MinMyGID();
2443 end = matrix->RowMap().MaxMyGID()+1;
2445 begin = matrix->RowMap().MinMyGID64();
2446 end = matrix->RowMap().MaxMyGID64()+1;
2449 return ((index >= static_cast<size_type>(begin)) &&
2450 (index < static_cast<size_type>(end)));
2471 const TrilinosScalar value)
2476 set (i, 1, &j, &value,
false);
2485 const bool elide_zero_values)
2487 Assert (indices.size() == values.
m(),
2488 ExcDimensionMismatch(indices.size(), values.
m()));
2489 Assert (values.
m() == values.
n(), ExcNotQuadratic());
2491 for (size_type i=0; i<indices.size(); ++i)
2492 set (indices[i], indices.size(), &indices[0], &values(i,0),
2502 const TrilinosScalar value)
2515 if (last_action == Insert)
2518 ierr = matrix->GlobalAssemble(*column_space_map,
2519 matrix->RowMap(),
false);
2521 Assert (ierr == 0, ExcTrilinosError(ierr));
2530 add (i, 1, &j, &value,
false);
2541 #ifndef DEAL_II_WITH_64BIT_INDICES 2542 return matrix->NumGlobalRows();
2544 return matrix->NumGlobalRows64();
2557 Assert(column_space_map.get() != 0, ExcInternalError());
2558 #ifndef DEAL_II_WITH_64BIT_INDICES 2559 return column_space_map->NumGlobalElements();
2561 return column_space_map->NumGlobalElements64();
2571 return matrix -> NumMyRows();
2577 std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
2580 size_type begin, end;
2581 #ifndef DEAL_II_WITH_64BIT_INDICES 2582 begin = matrix->RowMap().MinMyGID();
2583 end = matrix->RowMap().MaxMyGID()+1;
2585 begin = matrix->RowMap().MinMyGID64();
2586 end = matrix->RowMap().MaxMyGID64()+1;
2589 return std::make_pair (begin, end);
2598 #ifndef DEAL_II_WITH_64BIT_INDICES 2599 return matrix->NumGlobalNonzeros();
2601 return matrix->NumGlobalNonzeros64();
2607 template <
typename SparsityPatternType>
2610 const SparsityPatternType &sparsity_pattern,
2611 const MPI_Comm &communicator,
2612 const bool exchange_data)
2614 reinit (parallel_partitioning, parallel_partitioning,
2615 sparsity_pattern, communicator, exchange_data);
2620 template <
typename number>
2623 const ::SparseMatrix<number> &sparse_matrix,
2624 const MPI_Comm &communicator,
2625 const double drop_tolerance,
2626 const bool copy_values,
2627 const ::SparsityPattern *use_this_sparsity)
2630 reinit (parallel_partitioning, parallel_partitioning, sparse_matrix,
2631 drop_tolerance, copy_values, use_this_sparsity);
2637 const Epetra_CrsMatrix &
2640 return static_cast<const Epetra_CrsMatrix &
>(*matrix);
2646 const Epetra_CrsGraph &
2649 return matrix->Graph();
2658 return IndexSet(matrix->DomainMap());
2667 return IndexSet(matrix->RangeMap());
2695 DEAL_II_NAMESPACE_CLOSE
2698 #endif // DEAL_II_WITH_TRILINOS std_cxx11::shared_ptr< std::vector< TrilinosScalar > > value_cache
bool operator<(const Iterator< Constness > &) const
const_iterator end() const
bool operator>(const Iterator< Constness > &) const
std_cxx11::shared_ptr< std::vector< size_type > > colnum_cache
#define DeclException2(Exception2, type1, type2, outsequence)
SparseMatrixIterators::Iterator< false > iterator
::types::global_dof_index size_type
bool operator==(const Iterator< Constness > &) const
Accessor(MatrixType *matrix, const size_type row, const size_type index)
bool is_compressed() const
::types::global_dof_index size_type
DeclException0(ExcBeyondEndOfMatrix)
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
BlockDynamicSparsityPattern BlockCompressedSparsityPattern DEAL_II_DEPRECATED
bool in_local_range(const size_type index) const
const Accessor< Constness > & operator*() const
#define DeclException1(Exception1, type1, outsequence)
void add(const size_type i, const size_type j, const TrilinosScalar value)
#define Assert(cond, exc)
bool operator!=(const Iterator< Constness > &) const
const Accessor< Constness > * operator->() const
const Epetra_CrsMatrix & trilinos_matrix() const
AccessorBase(SparseMatrix *matrix, const size_type row, const size_type index)
const SparseMatrix MatrixType
IndexSet locally_owned_range_indices() const
size_type n_nonzero_elements() const
const_iterator begin() const
Accessor< Constness > accessor
SparseMatrixIterators::Iterator< true > const_iterator
void reinit(const SparsityPatternType &sparsity_pattern)
DeclException3(ExcAccessToNonlocalRow, std::size_t, std::size_t, std::size_t,<< "You tried to access row "<< arg1<< " of a distributed sparsity pattern, "<< " but only rows "<< arg2<< " through "<< arg3<< " are stored locally and can be accessed.")
Iterator< Constness > & operator++()
::types::global_dof_index size_type
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Accessor< Constness >::MatrixType MatrixType
std::pair< size_type, size_type > local_range() const
Accessor(MatrixType *matrix, const size_type row, const size_type index)
IndexSet locally_owned_domain_indices() const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
TrilinosScalar value() const
TrilinosScalar value_type
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
unsigned int local_size() const