Reference documentation for deal.II version 8.4.2
trilinos_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__trilinos_sparse_matrix_h
17 #define dealii__trilinos_sparse_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
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>
31 
32 # include <vector>
33 # include <cmath>
34 # include <memory>
35 
36 # define TrilinosScalar double
37 
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>
45 # include "mpi.h"
46 # else
47 # include "Epetra_SerialComm.h"
48 # endif
50 
51 class Epetra_Export;
52 
53 DEAL_II_NAMESPACE_OPEN
54 
55 // forward declarations
56 template <typename MatrixType> class BlockMatrixBase;
57 
58 template <typename number> class SparseMatrix;
59 class SparsityPattern;
61 
62 
63 namespace TrilinosWrappers
64 {
65  // forward declarations
66  class SparseMatrix;
67  class SparsityPattern;
68 
73  {
74  // forward declaration
75  template <bool Constness> class Iterator;
76 
80  DeclException0 (ExcBeyondEndOfMatrix);
81 
85  DeclException3 (ExcAccessToNonlocalRow,
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.");
91 
103  {
104  public:
108  typedef ::types::global_dof_index size_type;
109 
114  const size_type row,
115  const size_type index);
116 
120  size_type row() const;
121 
125  size_type index() const;
126 
130  size_type column() const;
131 
132  protected:
143  size_type a_row;
144 
148  size_type a_index;
149 
155  void visit_present_row ();
156 
169  std_cxx11::shared_ptr<std::vector<size_type> > colnum_cache;
170 
174  std_cxx11::shared_ptr<std::vector<TrilinosScalar> > value_cache;
175  };
176 
187  template <bool Constess>
188  class Accessor : public AccessorBase
189  {
193  TrilinosScalar value() const;
194 
198  TrilinosScalar &value();
199  };
200 
204  template<>
205  class Accessor<true> : public AccessorBase
206  {
207  public:
212  typedef const SparseMatrix MatrixType;
213 
218  Accessor (MatrixType *matrix,
219  const size_type row,
220  const size_type index);
221 
226  template <bool Other>
227  Accessor (const Accessor<Other> &a);
228 
232  TrilinosScalar value() const;
233 
234  private:
238  template <bool> friend class Iterator;
239  };
240 
244  template<>
245  class Accessor<false> : public AccessorBase
246  {
247  class Reference
248  {
249  public:
253  Reference (const Accessor<false> &accessor);
254 
258  operator TrilinosScalar () const;
259 
263  const Reference &operator = (const TrilinosScalar n) const;
264 
268  const Reference &operator += (const TrilinosScalar n) const;
269 
273  const Reference &operator -= (const TrilinosScalar n) const;
274 
278  const Reference &operator *= (const TrilinosScalar n) const;
279 
283  const Reference &operator /= (const TrilinosScalar n) const;
284 
285  private:
290  Accessor &accessor;
291  };
292 
293  public:
299 
304  Accessor (MatrixType *matrix,
305  const size_type row,
306  const size_type index);
307 
311  Reference value() const;
312 
313  private:
317  template <bool> friend class Iterator;
321  friend class Reference;
322  };
323 
338  template <bool Constness>
339  class Iterator
340  {
341  public:
345  typedef ::types::global_dof_index size_type;
346 
352 
357  Iterator (MatrixType *matrix,
358  const size_type row,
359  const size_type index);
360 
364  template <bool Other>
365  Iterator(const Iterator<Other> &other);
366 
370  Iterator<Constness> &operator++ ();
371 
375  Iterator<Constness> operator++ (int);
376 
380  const Accessor<Constness> &operator* () const;
381 
385  const Accessor<Constness> *operator-> () const;
386 
391  bool operator == (const Iterator<Constness> &) const;
392 
396  bool operator != (const Iterator<Constness> &) const;
397 
403  bool operator < (const Iterator<Constness> &) const;
404 
408  bool operator > (const Iterator<Constness> &) const;
409 
413  DeclException2 (ExcInvalidIndexWithinRow,
414  size_type, size_type,
415  << "Attempt to access element " << arg2
416  << " of row " << arg1
417  << " which doesn't have that many elements.");
418 
419  private:
424 
425  template <bool Other> friend class Iterator;
426  };
427 
428  }
429 
430 
492  class SparseMatrix : public Subscriptor
493  {
494  public:
498  typedef ::types::global_dof_index size_type;
499 
507  struct Traits
508  {
513  static const bool zero_addition_can_be_elided = true;
514  };
515 
520 
525 
529  typedef TrilinosScalar value_type;
530 
538  SparseMatrix ();
539 
547  SparseMatrix (const size_type m,
548  const size_type n,
549  const unsigned int n_max_entries_per_row);
550 
558  SparseMatrix (const size_type m,
559  const size_type n,
560  const std::vector<unsigned int> &n_entries_per_row);
561 
565  SparseMatrix (const SparsityPattern &InputSparsityPattern);
566 
570  virtual ~SparseMatrix ();
571 
587  template<typename SparsityPatternType>
588  void reinit (const SparsityPatternType &sparsity_pattern);
589 
602  void reinit (const SparsityPattern &sparsity_pattern);
603 
612  void reinit (const SparseMatrix &sparse_matrix);
613 
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);
639 
645  void reinit (const Epetra_CrsMatrix &input_matrix,
646  const bool copy_values = true);
648 
666  SparseMatrix (const Epetra_Map &parallel_partitioning,
667  const size_type n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
668 
678  SparseMatrix (const Epetra_Map &parallel_partitioning,
679  const std::vector<unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
680 
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;
702 
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;
722 
749  template<typename SparsityPatternType>
750  void reinit (const Epetra_Map &parallel_partitioning,
751  const SparsityPatternType &sparsity_pattern,
752  const bool exchange_data = false) DEAL_II_DEPRECATED;
753 
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;
773 
792  template <typename number>
793  void reinit (const Epetra_Map &parallel_partitioning,
794  const ::SparseMatrix<number> &dealii_sparse_matrix,
795  const double drop_tolerance=1e-13,
796  const bool copy_values=true,
797  const ::SparsityPattern *use_this_sparsity=0) DEAL_II_DEPRECATED;
798 
814  template <typename number>
815  void reinit (const Epetra_Map &row_parallel_partitioning,
816  const Epetra_Map &col_parallel_partitioning,
817  const ::SparseMatrix<number> &dealii_sparse_matrix,
818  const double drop_tolerance=1e-13,
819  const bool copy_values=true,
820  const ::SparsityPattern *use_this_sparsity=0) DEAL_II_DEPRECATED;
822 
838  SparseMatrix (const IndexSet &parallel_partitioning,
839  const MPI_Comm &communicator = MPI_COMM_WORLD,
840  const unsigned int n_max_entries_per_row = 0);
841 
849  SparseMatrix (const IndexSet &parallel_partitioning,
850  const MPI_Comm &communicator,
851  const std::vector<unsigned int> &n_entries_per_row);
852 
867  SparseMatrix (const IndexSet &row_parallel_partitioning,
868  const IndexSet &col_parallel_partitioning,
869  const MPI_Comm &communicator = MPI_COMM_WORLD,
870  const size_type n_max_entries_per_row = 0);
871 
886  SparseMatrix (const IndexSet &row_parallel_partitioning,
887  const IndexSet &col_parallel_partitioning,
888  const MPI_Comm &communicator,
889  const std::vector<unsigned int> &n_entries_per_row);
890 
911  template<typename SparsityPatternType>
912  void reinit (const IndexSet &parallel_partitioning,
913  const SparsityPatternType &sparsity_pattern,
914  const MPI_Comm &communicator = MPI_COMM_WORLD,
915  const bool exchange_data = false);
916 
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);
935 
952  template <typename number>
953  void reinit (const IndexSet &parallel_partitioning,
954  const ::SparseMatrix<number> &dealii_sparse_matrix,
955  const MPI_Comm &communicator = MPI_COMM_WORLD,
956  const double drop_tolerance=1e-13,
957  const bool copy_values=true,
958  const ::SparsityPattern *use_this_sparsity=0);
959 
973  template <typename number>
974  void reinit (const IndexSet &row_parallel_partitioning,
975  const IndexSet &col_parallel_partitioning,
976  const ::SparseMatrix<number> &dealii_sparse_matrix,
977  const MPI_Comm &communicator = MPI_COMM_WORLD,
978  const double drop_tolerance=1e-13,
979  const bool copy_values=true,
980  const ::SparsityPattern *use_this_sparsity=0);
982 
986 
990  size_type m () const;
991 
995  size_type n () const;
996 
1005  unsigned int local_size () const;
1006 
1015  std::pair<size_type, size_type>
1016  local_range () const;
1017 
1022  bool in_local_range (const size_type index) const;
1023 
1027  size_type n_nonzero_elements () const;
1028 
1032  unsigned int row_length (const size_type row) const;
1033 
1040  bool is_compressed () const;
1041 
1047  size_type memory_consumption () const;
1048 
1052  MPI_Comm get_mpi_communicator () const;
1053 
1055 
1059 
1069  SparseMatrix &
1070  operator = (const double d);
1071 
1079  void clear ();
1080 
1108  void compress (::VectorOperation::values operation);
1109 
1131  void set (const size_type i,
1132  const size_type j,
1133  const TrilinosScalar value);
1134 
1167  void set (const std::vector<size_type> &indices,
1168  const FullMatrix<TrilinosScalar> &full_matrix,
1169  const bool elide_zero_values = false);
1170 
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);
1180 
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);
1212 
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);
1245 
1255  void add (const size_type i,
1256  const size_type j,
1257  const TrilinosScalar value);
1258 
1277  void add (const std::vector<size_type> &indices,
1278  const FullMatrix<TrilinosScalar> &full_matrix,
1279  const bool elide_zero_values = true);
1280 
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);
1290 
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);
1308 
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);
1328 
1332  SparseMatrix &operator *= (const TrilinosScalar factor);
1333 
1337  SparseMatrix &operator /= (const TrilinosScalar factor);
1338 
1342  void copy_from (const SparseMatrix &source);
1343 
1351  void add (const TrilinosScalar factor,
1352  const SparseMatrix &matrix);
1353 
1380  void clear_row (const size_type row,
1381  const TrilinosScalar new_diag_value = 0);
1382 
1403  void clear_rows (const std::vector<size_type> &rows,
1404  const TrilinosScalar new_diag_value = 0);
1405 
1412  void transpose ();
1413 
1415 
1419 
1428  TrilinosScalar operator () (const size_type i,
1429  const size_type j) const;
1430 
1447  TrilinosScalar el (const size_type i,
1448  const size_type j) const;
1449 
1456  TrilinosScalar diag_element (const size_type i) const;
1457 
1459 
1463 
1485  template<typename VectorType>
1486  void vmult (VectorType &dst,
1487  const VectorType &src) const;
1488 
1511  template <typename VectorType>
1512  void Tvmult (VectorType &dst,
1513  const VectorType &src) const;
1514 
1538  template<typename VectorType>
1539  void vmult_add (VectorType &dst,
1540  const VectorType &src) const;
1541 
1565  template <typename VectorType>
1566  void Tvmult_add (VectorType &dst,
1567  const VectorType &src) const;
1568 
1593  TrilinosScalar matrix_norm_square (const VectorBase &v) const;
1594 
1613  TrilinosScalar matrix_scalar_product (const VectorBase &u,
1614  const VectorBase &v) const;
1615 
1633  TrilinosScalar residual (VectorBase &dst,
1634  const VectorBase &x,
1635  const VectorBase &b) const;
1636 
1651  void mmult (SparseMatrix &C,
1652  const SparseMatrix &B,
1653  const VectorBase &V = VectorBase()) const;
1654 
1655 
1672  void Tmmult (SparseMatrix &C,
1673  const SparseMatrix &B,
1674  const VectorBase &V = VectorBase()) const;
1675 
1677 
1681 
1689  TrilinosScalar l1_norm () const;
1690 
1699  TrilinosScalar linfty_norm () const;
1700 
1705  TrilinosScalar frobenius_norm () const;
1706 
1708 
1712 
1717  const Epetra_CrsMatrix &trilinos_matrix () const;
1718 
1723  const Epetra_CrsGraph &trilinos_sparsity_pattern () const;
1724 
1732  const Epetra_Map &domain_partitioner () const DEAL_II_DEPRECATED;
1733 
1742  const Epetra_Map &range_partitioner () const DEAL_II_DEPRECATED;
1743 
1751  const Epetra_Map &row_partitioner () const DEAL_II_DEPRECATED;
1752 
1762  const Epetra_Map &col_partitioner () const DEAL_II_DEPRECATED;
1764 
1769 
1774  IndexSet locally_owned_domain_indices() const;
1775 
1781  IndexSet locally_owned_range_indices() const;
1782 
1784 
1789 
1808  const_iterator begin () const;
1809 
1813  iterator begin ();
1814 
1819  const_iterator end () const;
1820 
1824  iterator end ();
1825 
1854  const_iterator begin (const size_type r) const;
1855 
1859  iterator begin (const size_type r);
1860 
1870  const_iterator end (const size_type r) const;
1871 
1875  iterator end (const size_type r);
1876 
1878 
1882 
1888  void write_ascii ();
1889 
1897  void print (std::ostream &out,
1898  const bool write_extended_trilinos_info = false) const;
1899 
1901 
1909  DeclException1 (ExcTrilinosError,
1910  int,
1911  << "An error with error number " << arg1
1912  << " occurred while calling a Trilinos function");
1913 
1917  DeclException2 (ExcInvalidIndex,
1918  size_type, size_type,
1919  << "The entry with index <" << arg1 << ',' << arg2
1920  << "> does not exist.");
1921 
1925  DeclException0 (ExcSourceEqualsDestination);
1926 
1930  DeclException0 (ExcMatrixNotCompressed);
1931 
1935  DeclException4 (ExcAccessToNonLocalElement,
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.");
1942 
1946  DeclException2 (ExcAccessToNonPresentElement,
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.");
1953 
1954 
1955 
1956  protected:
1957 
1968  void prepare_add();
1969 
1975  void prepare_set();
1976 
1977 
1978 
1979  private:
1983  SparseMatrix (const SparseMatrix &);
1987  SparseMatrix &operator = (const SparseMatrix &);
1988 
1993  std_cxx11::shared_ptr<Epetra_Map> column_space_map;
1994 
2000  std_cxx11::shared_ptr<Epetra_FECrsMatrix> matrix;
2001 
2007  std_cxx11::shared_ptr<Epetra_CrsMatrix> nonlocal_matrix;
2008 
2012  std_cxx11::shared_ptr<Epetra_Export> nonlocal_matrix_exporter;
2013 
2025  Epetra_CombineMode last_action;
2026 
2031  bool compressed;
2032 
2036  friend class BlockMatrixBase<SparseMatrix>;
2037  };
2038 
2039 
2040 
2041 // -------------------------- inline and template functions ----------------------
2042 
2043 
2044 #ifndef DOXYGEN
2045 
2046  namespace SparseMatrixIterators
2047  {
2048  inline
2049  AccessorBase::AccessorBase(SparseMatrix *matrix, size_type row, size_type index)
2050  :
2051  matrix(matrix),
2052  a_row(row),
2053  a_index(index)
2054  {
2055  visit_present_row ();
2056  }
2057 
2058 
2059  inline
2061  AccessorBase::row() const
2062  {
2063  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2064  return a_row;
2065  }
2066 
2067 
2068  inline
2070  AccessorBase::column() const
2071  {
2072  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2073  return (*colnum_cache)[a_index];
2074  }
2075 
2076 
2077  inline
2079  AccessorBase::index() const
2080  {
2081  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2082  return a_index;
2083  }
2084 
2085 
2086  inline
2087  Accessor<true>::Accessor (MatrixType *matrix,
2088  const size_type row,
2089  const size_type index)
2090  :
2091  AccessorBase(const_cast<SparseMatrix *>(matrix), row, index)
2092  {}
2093 
2094 
2095  template <bool Other>
2096  inline
2097  Accessor<true>::Accessor(const Accessor<Other> &other)
2098  :
2099  AccessorBase(other)
2100  {}
2101 
2102 
2103  inline
2104  TrilinosScalar
2105  Accessor<true>::value() const
2106  {
2107  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2108  return (*value_cache)[a_index];
2109  }
2110 
2111 
2112  inline
2113  Accessor<false>::Reference::Reference (
2114  const Accessor<false> &acc)
2115  :
2116  accessor(const_cast<Accessor<false>&>(acc))
2117  {}
2118 
2119 
2120  inline
2121  Accessor<false>::Reference::operator TrilinosScalar () const
2122  {
2123  return (*accessor.value_cache)[accessor.a_index];
2124  }
2125 
2126  inline
2128  Accessor<false>::Reference::operator = (const TrilinosScalar n) const
2129  {
2130  (*accessor.value_cache)[accessor.a_index] = n;
2131  accessor.matrix->set(accessor.row(), accessor.column(),
2132  static_cast<TrilinosScalar>(*this));
2133  return *this;
2134  }
2135 
2136 
2137  inline
2139  Accessor<false>::Reference::operator += (const TrilinosScalar n) const
2140  {
2141  (*accessor.value_cache)[accessor.a_index] += n;
2142  accessor.matrix->set(accessor.row(), accessor.column(),
2143  static_cast<TrilinosScalar>(*this));
2144  return *this;
2145  }
2146 
2147 
2148  inline
2150  Accessor<false>::Reference::operator -= (const TrilinosScalar n) const
2151  {
2152  (*accessor.value_cache)[accessor.a_index] -= n;
2153  accessor.matrix->set(accessor.row(), accessor.column(),
2154  static_cast<TrilinosScalar>(*this));
2155  return *this;
2156  }
2157 
2158 
2159  inline
2161  Accessor<false>::Reference::operator *= (const TrilinosScalar n) const
2162  {
2163  (*accessor.value_cache)[accessor.a_index] *= n;
2164  accessor.matrix->set(accessor.row(), accessor.column(),
2165  static_cast<TrilinosScalar>(*this));
2166  return *this;
2167  }
2168 
2169 
2170  inline
2172  Accessor<false>::Reference::operator /= (const TrilinosScalar n) const
2173  {
2174  (*accessor.value_cache)[accessor.a_index] /= n;
2175  accessor.matrix->set(accessor.row(), accessor.column(),
2176  static_cast<TrilinosScalar>(*this));
2177  return *this;
2178  }
2179 
2180 
2181  inline
2182  Accessor<false>::Accessor (MatrixType *matrix,
2183  const size_type row,
2184  const size_type index)
2185  :
2186  AccessorBase(matrix, row, index)
2187  {}
2188 
2189 
2190  inline
2192  Accessor<false>::value() const
2193  {
2194  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2195  return Reference(*this);
2196  }
2197 
2198 
2199 
2200  template <bool Constness>
2201  inline
2202  Iterator<Constness>::Iterator(MatrixType *matrix,
2203  const size_type row,
2204  const size_type index)
2205  :
2206  accessor(matrix, row, index)
2207  {}
2208 
2209 
2210  template <bool Constness>
2211  template <bool Other>
2212  inline
2213  Iterator<Constness>::Iterator(const Iterator<Other> &other)
2214  :
2215  accessor(other.accessor)
2216  {}
2217 
2218 
2219  template <bool Constness>
2220  inline
2221  Iterator<Constness> &
2223  {
2224  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
2225 
2226  ++accessor.a_index;
2227 
2228  // If at end of line: do one
2229  // step, then cycle until we
2230  // find a row with a nonzero
2231  // number of entries.
2232  if (accessor.a_index >= accessor.colnum_cache->size())
2233  {
2234  accessor.a_index = 0;
2235  ++accessor.a_row;
2236 
2237  while ((accessor.a_row < accessor.matrix->m())
2238  &&
2239  ((accessor.matrix->in_local_range (accessor.a_row) == false)
2240  ||
2241  (accessor.matrix->row_length(accessor.a_row) == 0)))
2242  ++accessor.a_row;
2243 
2244  accessor.visit_present_row();
2245  }
2246  return *this;
2247  }
2248 
2249 
2250  template <bool Constness>
2251  inline
2252  Iterator<Constness>
2254  {
2255  const Iterator<Constness> old_state = *this;
2256  ++(*this);
2257  return old_state;
2258  }
2259 
2260 
2261 
2262  template <bool Constness>
2263  inline
2264  const Accessor<Constness> &
2266  {
2267  return accessor;
2268  }
2269 
2270 
2271 
2272  template <bool Constness>
2273  inline
2274  const Accessor<Constness> *
2276  {
2277  return &accessor;
2278  }
2279 
2280 
2281 
2282  template <bool Constness>
2283  inline
2284  bool
2285  Iterator<Constness>::operator == (const Iterator<Constness> &other) const
2286  {
2287  return (accessor.a_row == other.accessor.a_row &&
2288  accessor.a_index == other.accessor.a_index);
2289  }
2290 
2291 
2292 
2293  template <bool Constness>
2294  inline
2295  bool
2296  Iterator<Constness>::operator != (const Iterator<Constness> &other) const
2297  {
2298  return ! (*this == other);
2299  }
2300 
2301 
2302 
2303  template <bool Constness>
2304  inline
2305  bool
2306  Iterator<Constness>::operator < (const Iterator<Constness> &other) const
2307  {
2308  return (accessor.row() < other.accessor.row() ||
2309  (accessor.row() == other.accessor.row() &&
2310  accessor.index() < other.accessor.index()));
2311  }
2312 
2313 
2314  template <bool Constness>
2315  inline
2316  bool
2317  Iterator<Constness>::operator > (const Iterator<Constness> &other) const
2318  {
2319  return (other < *this);
2320  }
2321 
2322  }
2323 
2324 
2325 
2326  inline
2328  SparseMatrix::begin() const
2329  {
2330  return begin(0);
2331  }
2332 
2333 
2334 
2335  inline
2337  SparseMatrix::end() const
2338  {
2339  return const_iterator(this, m(), 0);
2340  }
2341 
2342 
2343 
2344  inline
2346  SparseMatrix::begin(const size_type r) const
2347  {
2348  Assert (r < m(), ExcIndexRange(r, 0, m()));
2349  if (in_local_range (r)
2350  &&
2351  (row_length(r) > 0))
2352  return const_iterator(this, r, 0);
2353  else
2354  return end (r);
2355  }
2356 
2357 
2358 
2359  inline
2361  SparseMatrix::end(const size_type r) const
2362  {
2363  Assert (r < m(), ExcIndexRange(r, 0, m()));
2364 
2365  // place the iterator on the first entry
2366  // past this line, or at the end of the
2367  // matrix
2368  for (size_type i=r+1; i<m(); ++i)
2369  if (in_local_range (i)
2370  &&
2371  (row_length(i) > 0))
2372  return const_iterator(this, i, 0);
2373 
2374  // if there is no such line, then take the
2375  // end iterator of the matrix
2376  return end();
2377  }
2378 
2379 
2380 
2381  inline
2384  {
2385  return begin(0);
2386  }
2387 
2388 
2389 
2390  inline
2393  {
2394  return iterator(this, m(), 0);
2395  }
2396 
2397 
2398 
2399  inline
2401  SparseMatrix::begin(const size_type r)
2402  {
2403  Assert (r < m(), ExcIndexRange(r, 0, m()));
2404  if (in_local_range (r)
2405  &&
2406  (row_length(r) > 0))
2407  return iterator(this, r, 0);
2408  else
2409  return end (r);
2410  }
2411 
2412 
2413 
2414  inline
2416  SparseMatrix::end(const size_type r)
2417  {
2418  Assert (r < m(), ExcIndexRange(r, 0, m()));
2419 
2420  // place the iterator on the first entry
2421  // past this line, or at the end of the
2422  // matrix
2423  for (size_type i=r+1; i<m(); ++i)
2424  if (in_local_range (i)
2425  &&
2426  (row_length(i) > 0))
2427  return iterator(this, i, 0);
2428 
2429  // if there is no such line, then take the
2430  // end iterator of the matrix
2431  return end();
2432  }
2433 
2434 
2435 
2436  inline
2437  bool
2438  SparseMatrix::in_local_range (const size_type index) const
2439  {
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;
2444 #else
2445  begin = matrix->RowMap().MinMyGID64();
2446  end = matrix->RowMap().MaxMyGID64()+1;
2447 #endif
2448 
2449  return ((index >= static_cast<size_type>(begin)) &&
2450  (index < static_cast<size_type>(end)));
2451  }
2452 
2453 
2454 
2455  inline
2456  bool
2458  {
2459  return compressed;
2460  }
2461 
2462 
2463 
2464  // Inline the set() and add() functions, since they will be called
2465  // frequently, and the compiler can optimize away some unnecessary loops
2466  // when the sizes are given at compile time.
2467  inline
2468  void
2469  SparseMatrix::set (const size_type i,
2470  const size_type j,
2471  const TrilinosScalar value)
2472  {
2473 
2474  AssertIsFinite(value);
2475 
2476  set (i, 1, &j, &value, false);
2477  }
2478 
2479 
2480 
2481  inline
2482  void
2483  SparseMatrix::set (const std::vector<size_type> &indices,
2484  const FullMatrix<TrilinosScalar> &values,
2485  const bool elide_zero_values)
2486  {
2487  Assert (indices.size() == values.m(),
2488  ExcDimensionMismatch(indices.size(), values.m()));
2489  Assert (values.m() == values.n(), ExcNotQuadratic());
2490 
2491  for (size_type i=0; i<indices.size(); ++i)
2492  set (indices[i], indices.size(), &indices[0], &values(i,0),
2493  elide_zero_values);
2494  }
2495 
2496 
2497 
2498  inline
2499  void
2500  SparseMatrix::add (const size_type i,
2501  const size_type j,
2502  const TrilinosScalar value)
2503  {
2504  AssertIsFinite(value);
2505 
2506  if (value == 0)
2507  {
2508  // we have to check after Insert/Add in any case to be consistent
2509  // with the MPI communication model (see the comments in the
2510  // documentation of TrilinosWrappers::Vector), but we can save some
2511  // work if the addend is zero. However, these actions are done in case
2512  // we pass on to the other function.
2513 
2514  // TODO: fix this (do not run compress here, but fail)
2515  if (last_action == Insert)
2516  {
2517  int ierr;
2518  ierr = matrix->GlobalAssemble(*column_space_map,
2519  matrix->RowMap(), false);
2520 
2521  Assert (ierr == 0, ExcTrilinosError(ierr));
2522  (void)ierr; // removes -Wunused-but-set-variable in optimized mode
2523  }
2524 
2525  last_action = Add;
2526 
2527  return;
2528  }
2529  else
2530  add (i, 1, &j, &value, false);
2531  }
2532 
2533 
2534 
2535  // inline "simple" functions that are called frequently and do only involve
2536  // a call to some Trilinos function.
2537  inline
2539  SparseMatrix::m () const
2540  {
2541 #ifndef DEAL_II_WITH_64BIT_INDICES
2542  return matrix->NumGlobalRows();
2543 #else
2544  return matrix->NumGlobalRows64();
2545 #endif
2546  }
2547 
2548 
2549 
2550  inline
2552  SparseMatrix::n () const
2553  {
2554  // If the matrix structure has not been fixed (i.e., we did not have a
2555  // sparsity pattern), it does not know about the number of columns so we
2556  // must always take this from the additional column space map
2557  Assert(column_space_map.get() != 0, ExcInternalError());
2558 #ifndef DEAL_II_WITH_64BIT_INDICES
2559  return column_space_map->NumGlobalElements();
2560 #else
2561  return column_space_map->NumGlobalElements64();
2562 #endif
2563  }
2564 
2565 
2566 
2567  inline
2568  unsigned int
2569  SparseMatrix::local_size () const
2570  {
2571  return matrix -> NumMyRows();
2572  }
2573 
2574 
2575 
2576  inline
2577  std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
2578  SparseMatrix::local_range () const
2579  {
2580  size_type begin, end;
2581 #ifndef DEAL_II_WITH_64BIT_INDICES
2582  begin = matrix->RowMap().MinMyGID();
2583  end = matrix->RowMap().MaxMyGID()+1;
2584 #else
2585  begin = matrix->RowMap().MinMyGID64();
2586  end = matrix->RowMap().MaxMyGID64()+1;
2587 #endif
2588 
2589  return std::make_pair (begin, end);
2590  }
2591 
2592 
2593 
2594  inline
2597  {
2598 #ifndef DEAL_II_WITH_64BIT_INDICES
2599  return matrix->NumGlobalNonzeros();
2600 #else
2601  return matrix->NumGlobalNonzeros64();
2602 #endif
2603  }
2604 
2605 
2606 
2607  template <typename SparsityPatternType>
2608  inline
2609  void SparseMatrix::reinit (const IndexSet &parallel_partitioning,
2610  const SparsityPatternType &sparsity_pattern,
2611  const MPI_Comm &communicator,
2612  const bool exchange_data)
2613  {
2614  reinit (parallel_partitioning, parallel_partitioning,
2615  sparsity_pattern, communicator, exchange_data);
2616  }
2617 
2618 
2619 
2620  template <typename number>
2621  inline
2622  void SparseMatrix::reinit (const IndexSet &parallel_partitioning,
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)
2628  {
2629  Epetra_Map map = parallel_partitioning.make_trilinos_map (communicator, false);
2630  reinit (parallel_partitioning, parallel_partitioning, sparse_matrix,
2631  drop_tolerance, copy_values, use_this_sparsity);
2632  }
2633 
2634 
2635 
2636  inline
2637  const Epetra_CrsMatrix &
2639  {
2640  return static_cast<const Epetra_CrsMatrix &>(*matrix);
2641  }
2642 
2643 
2644 
2645  inline
2646  const Epetra_CrsGraph &
2648  {
2649  return matrix->Graph();
2650  }
2651 
2652 
2653 
2654  inline
2655  IndexSet
2657  {
2658  return IndexSet(matrix->DomainMap());
2659  }
2660 
2661 
2662 
2663  inline
2664  IndexSet
2666  {
2667  return IndexSet(matrix->RangeMap());
2668  }
2669 
2670 
2671 
2672  inline
2673  void
2675  {
2676  //nothing to do here
2677  }
2678 
2679 
2680 
2681  inline
2682  void
2684  {
2685  //nothing to do here
2686  }
2687 
2688 
2689 
2690 #endif // DOXYGEN
2691 
2692 }
2693 
2694 
2695 DEAL_II_NAMESPACE_CLOSE
2696 
2697 
2698 #endif // DEAL_II_WITH_TRILINOS
2699 
2700 
2701 /*----------------------- trilinos_sparse_matrix.h --------------------*/
2702 
2703 #endif
2704 /*----------------------- trilinos_sparse_matrix.h --------------------*/
std_cxx11::shared_ptr< std::vector< TrilinosScalar > > value_cache
bool operator<(const Iterator< Constness > &) const
size_type m() 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)
Definition: exceptions.h:552
SparseMatrixIterators::Iterator< false > iterator
bool operator==(const Iterator< Constness > &) const
Accessor(MatrixType *matrix, const size_type row, const size_type index)
::types::global_dof_index size_type
STL namespace.
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:484
BlockDynamicSparsityPattern BlockCompressedSparsityPattern DEAL_II_DEPRECATED
bool in_local_range(const size_type index) const
const Accessor< Constness > & operator*() const
size_type n() const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
void add(const size_type i, const size_type j, const TrilinosScalar value)
#define Assert(cond, exc)
Definition: exceptions.h:294
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)
IndexSet locally_owned_range_indices() const
size_type n_nonzero_elements() const
const_iterator begin() const
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.")
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:572
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
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
Definition: exceptions.h:1096
unsigned int local_size() const