Reference documentation for deal.II version 8.5.1
sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2017 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__sparse_matrix_h
17 #define dealii__sparse_matrix_h
18 
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/subscriptor.h>
22 #include <deal.II/base/smartpointer.h>
23 #include <deal.II/lac/sparsity_pattern.h>
24 #include <deal.II/lac/identity_matrix.h>
25 #include <deal.II/lac/exceptions.h>
26 #include <deal.II/lac/vector.h>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 template <typename number> class Vector;
31 template <typename number> class FullMatrix;
32 template <typename Matrix> class BlockMatrixBase;
33 template <typename number> class SparseILU;
34 
35 #ifdef DEAL_II_WITH_TRILINOS
36 namespace TrilinosWrappers
37 {
38  class SparseMatrix;
39 }
40 #endif
41 
52 {
57 
58  // forward declaration
59  template <typename number, bool Constness>
60  class Iterator;
61 
72  template <typename number, bool Constness>
74  {
75  public:
79  number value() const;
80 
84  number &value();
85 
90  const SparseMatrix<number> &get_matrix () const;
91  };
92 
93 
94 
101  template <typename number>
102  class Accessor<number,true> : public SparsityPatternIterators::Accessor
103  {
104  public:
110 
114  Accessor (MatrixType *matrix,
115  const std::size_t index_within_matrix);
116 
120  Accessor (MatrixType *matrix);
121 
126 
130  number value() const;
131 
136  const MatrixType &get_matrix () const;
137 
138  private:
143 
148 
152  template <typename, bool>
153  friend class Iterator;
154  };
155 
156 
163  template <typename number>
164  class Accessor<number,false> : public SparsityPatternIterators::Accessor
165  {
166  private:
191  class Reference
192  {
193  public:
198  Reference (const Accessor *accessor,
199  const bool dummy);
200 
204  operator number () const;
205 
209  const Reference &operator = (const number n) const;
210 
214  const Reference &operator += (const number n) const;
215 
219  const Reference &operator -= (const number n) const;
220 
224  const Reference &operator *= (const number n) const;
225 
229  const Reference &operator /= (const number n) const;
230 
231  private:
237  };
238 
239  public:
245 
249  Accessor (MatrixType *matrix,
250  const std::size_t index);
251 
255  Accessor (MatrixType *matrix);
256 
260  Reference value() const;
261 
266  MatrixType &get_matrix () const;
267 
268  private:
273 
278 
282  template <typename, bool>
283  friend class Iterator;
284  };
285 
286 
287 
317  template <typename number, bool Constness>
318  class Iterator
319  {
320  public:
324  typedef
327 
332  typedef
334 
339  Iterator (MatrixType *matrix,
340  const std::size_t index_within_matrix);
341 
345  Iterator (MatrixType *matrix);
346 
352 
357 
361  Iterator operator++ (int);
362 
366  const Accessor<number,Constness> &operator* () const;
367 
372 
376  bool operator == (const Iterator &) const;
377 
381  bool operator != (const Iterator &) const;
382 
390  bool operator < (const Iterator &) const;
391 
396  bool operator > (const Iterator &) const;
397 
404  int operator - (const Iterator &p) const;
405 
409  Iterator operator + (const size_type n) const;
410 
411  private:
416  };
417 
418 }
419 
425 //TODO: Add multithreading to the other vmult functions.
426 
455 template <typename number>
456 class SparseMatrix : public virtual Subscriptor
457 {
458 public:
463 
468  typedef number value_type;
469 
480 
485  typedef
488 
495  typedef
498 
505  struct Traits
506  {
511  static const bool zero_addition_can_be_elided = true;
512  };
513 
528  SparseMatrix ();
529 
538  SparseMatrix (const SparseMatrix &);
539 
540 #ifdef DEAL_II_WITH_CXX11
541 
552 #endif
553 
567  explicit SparseMatrix (const SparsityPattern &sparsity);
568 
575  SparseMatrix (const SparsityPattern &sparsity,
576  const IdentityMatrix &id);
577 
582  virtual ~SparseMatrix ();
583 
594 
595 #ifdef DEAL_II_WITH_CXX11
596 
604 #endif
605 
613  operator= (const IdentityMatrix &id);
614 
626  SparseMatrix &operator = (const double d);
627 
641  virtual void reinit (const SparsityPattern &sparsity);
642 
648  virtual void clear ();
650 
658  bool empty () const;
659 
664  size_type m () const;
665 
670  size_type n () const;
671 
675  size_type get_row_length (const size_type row) const;
676 
682  std::size_t n_nonzero_elements () const;
683 
693  std::size_t n_actually_nonzero_elements (const double threshold = 0.) const;
694 
703  const SparsityPattern &get_sparsity_pattern () const;
704 
709  std::size_t memory_consumption () const;
710 
715 
717 
726  void set (const size_type i,
727  const size_type j,
728  const number value);
729 
745  template <typename number2>
746  void set (const std::vector<size_type> &indices,
747  const FullMatrix<number2> &full_matrix,
748  const bool elide_zero_values = false);
749 
755  template <typename number2>
756  void set (const std::vector<size_type> &row_indices,
757  const std::vector<size_type> &col_indices,
758  const FullMatrix<number2> &full_matrix,
759  const bool elide_zero_values = false);
760 
771  template <typename number2>
772  void set (const size_type row,
773  const std::vector<size_type> &col_indices,
774  const std::vector<number2> &values,
775  const bool elide_zero_values = false);
776 
786  template <typename number2>
787  void set (const size_type row,
788  const size_type n_cols,
789  const size_type *col_indices,
790  const number2 *values,
791  const bool elide_zero_values = false);
792 
798  void add (const size_type i,
799  const size_type j,
800  const number value);
801 
816  template <typename number2>
817  void add (const std::vector<size_type> &indices,
818  const FullMatrix<number2> &full_matrix,
819  const bool elide_zero_values = true);
820 
826  template <typename number2>
827  void add (const std::vector<size_type> &row_indices,
828  const std::vector<size_type> &col_indices,
829  const FullMatrix<number2> &full_matrix,
830  const bool elide_zero_values = true);
831 
841  template <typename number2>
842  void add (const size_type row,
843  const std::vector<size_type> &col_indices,
844  const std::vector<number2> &values,
845  const bool elide_zero_values = true);
846 
856  template <typename number2>
857  void add (const size_type row,
858  const size_type n_cols,
859  const size_type *col_indices,
860  const number2 *values,
861  const bool elide_zero_values = true,
862  const bool col_indices_are_sorted = false);
863 
867  SparseMatrix &operator *= (const number factor);
868 
872  SparseMatrix &operator /= (const number factor);
873 
886  void symmetrize ();
887 
904  template <typename somenumber>
906  copy_from (const SparseMatrix<somenumber> &source);
907 
924  template <typename ForwardIterator>
925  void copy_from (const ForwardIterator begin,
926  const ForwardIterator end);
927 
933  template <typename somenumber>
934  void copy_from (const FullMatrix<somenumber> &matrix);
935 
936 #ifdef DEAL_II_WITH_TRILINOS
937 
948 #endif
949 
961  template <typename somenumber>
962  void add (const number factor,
963  const SparseMatrix<somenumber> &matrix);
964 
966 
970 
984  number operator () (const size_type i,
985  const size_type j) const;
986 
999  number el (const size_type i,
1000  const size_type j) const;
1001 
1011  number diag_element (const size_type i) const;
1012 
1017  number &diag_element (const size_type i);
1018 
1020 
1040  template <class OutVector, class InVector>
1041  void vmult (OutVector &dst,
1042  const InVector &src) const;
1043 
1059  template <class OutVector, class InVector>
1060  void Tvmult (OutVector &dst,
1061  const InVector &src) const;
1062 
1079  template <class OutVector, class InVector>
1080  void vmult_add (OutVector &dst,
1081  const InVector &src) const;
1082 
1098  template <class OutVector, class InVector>
1099  void Tvmult_add (OutVector &dst,
1100  const InVector &src) const;
1101 
1119  template <typename somenumber>
1120  somenumber matrix_norm_square (const Vector<somenumber> &v) const;
1121 
1127  template <typename somenumber>
1128  somenumber matrix_scalar_product (const Vector<somenumber> &u,
1129  const Vector<somenumber> &v) const;
1130 
1140  template <typename somenumber>
1141  somenumber residual (Vector<somenumber> &dst,
1142  const Vector<somenumber> &x,
1143  const Vector<somenumber> &b) const;
1144 
1180  template <typename numberB, typename numberC>
1181  void mmult (SparseMatrix<numberC> &C,
1182  const SparseMatrix<numberB> &B,
1183  const Vector<number> &V = Vector<number>(),
1184  const bool rebuild_sparsity_pattern = true) const;
1185 
1210  template <typename numberB, typename numberC>
1211  void Tmmult (SparseMatrix<numberC> &C,
1212  const SparseMatrix<numberB> &B,
1213  const Vector<number> &V = Vector<number>(),
1214  const bool rebuild_sparsity_pattern = true) const;
1215 
1217 
1221 
1229  real_type l1_norm () const;
1230 
1238  real_type linfty_norm () const;
1239 
1244  real_type frobenius_norm () const;
1246 
1250 
1256  template <typename somenumber>
1258  const Vector<somenumber> &src,
1259  const number omega = 1.) const;
1260 
1267  template <typename somenumber>
1269  const Vector<somenumber> &src,
1270  const number omega = 1.,
1271  const std::vector<std::size_t> &pos_right_of_diagonal=std::vector<std::size_t>()) const;
1272 
1276  template <typename somenumber>
1278  const Vector<somenumber> &src,
1279  const number om = 1.) const;
1280 
1284  template <typename somenumber>
1286  const Vector<somenumber> &src,
1287  const number om = 1.) const;
1288 
1294  template <typename somenumber>
1295  void SSOR (Vector<somenumber> &v,
1296  const number omega = 1.) const;
1297 
1302  template <typename somenumber>
1303  void SOR (Vector<somenumber> &v,
1304  const number om = 1.) const;
1305 
1310  template <typename somenumber>
1311  void TSOR (Vector<somenumber> &v,
1312  const number om = 1.) const;
1313 
1324  template <typename somenumber>
1325  void PSOR (Vector<somenumber> &v,
1326  const std::vector<size_type> &permutation,
1327  const std::vector<size_type> &inverse_permutation,
1328  const number om = 1.) const;
1329 
1340  template <typename somenumber>
1341  void TPSOR (Vector<somenumber> &v,
1342  const std::vector<size_type> &permutation,
1343  const std::vector<size_type> &inverse_permutation,
1344  const number om = 1.) const;
1345 
1351  template <typename somenumber>
1353  const Vector<somenumber> &b,
1354  const number om = 1.) const;
1355 
1360  template <typename somenumber>
1361  void SOR_step (Vector<somenumber> &v,
1362  const Vector<somenumber> &b,
1363  const number om = 1.) const;
1364 
1369  template <typename somenumber>
1370  void TSOR_step (Vector<somenumber> &v,
1371  const Vector<somenumber> &b,
1372  const number om = 1.) const;
1373 
1378  template <typename somenumber>
1379  void SSOR_step (Vector<somenumber> &v,
1380  const Vector<somenumber> &b,
1381  const number om = 1.) const;
1383 
1387 
1394  const_iterator begin () const;
1395 
1399  iterator begin ();
1400 
1404  const_iterator end () const;
1405 
1409  iterator end ();
1410 
1420  const_iterator begin (const size_type r) const;
1421 
1425  iterator begin (const size_type r);
1426 
1436  const_iterator end (const size_type r) const;
1437 
1441  iterator end (const size_type r);
1443 
1447 
1459  template <class StreamType>
1460  void print (StreamType &out,
1461  const bool across = false,
1462  const bool diagonal_first = true) const;
1463 
1484  void print_formatted (std::ostream &out,
1485  const unsigned int precision = 3,
1486  const bool scientific = true,
1487  const unsigned int width = 0,
1488  const char *zero_string = " ",
1489  const double denominator = 1.) const;
1490 
1496  void print_pattern(std::ostream &out,
1497  const double threshold = 0.) const;
1498 
1509  void block_write (std::ostream &out) const;
1510 
1527  void block_read (std::istream &in);
1529 
1538  int, int,
1539  << "You are trying to access the matrix entry with index <"
1540  << arg1 << ',' << arg2
1541  << ">, but this entry does not exist in the sparsity pattern "
1542  "of this matrix."
1543  "\n\n"
1544  "The most common cause for this problem is that you used "
1545  "a method to build the sparsity pattern that did not "
1546  "(completely) take into account all of the entries you "
1547  "will later try to write into. An example would be "
1548  "building a sparsity pattern that does not include "
1549  "the entries you will write into due to constraints "
1550  "on degrees of freedom such as hanging nodes or periodic "
1551  "boundary conditions. In such cases, building the "
1552  "sparsity pattern will succeed, but you will get errors "
1553  "such as the current one at one point or other when "
1554  "trying to write into the entries of the matrix.");
1559  "When copying one sparse matrix into another, "
1560  "or when adding one sparse matrix to another, "
1561  "both matrices need to refer to the same "
1562  "sparsity pattern.");
1567  int, int,
1568  << "The iterators denote a range of " << arg1
1569  << " elements, but the given number of rows was " << arg2);
1575 
1576 protected:
1587  void prepare_add();
1588 
1593  void prepare_set();
1594 
1595 private:
1602 
1609  number *val;
1610 
1617  std::size_t max_len;
1618 
1619  // make all other sparse matrices friends
1620  template <typename somenumber> friend class SparseMatrix;
1621  template <typename somenumber> friend class SparseLUDecomposition;
1622  template <typename> friend class SparseILU;
1623 
1627  template <typename> friend class BlockMatrixBase;
1628 
1632  template <typename,bool> friend class SparseMatrixIterators::Iterator;
1633  template <typename,bool> friend class SparseMatrixIterators::Accessor;
1634 };
1635 
1636 #ifndef DOXYGEN
1637 /*---------------------- Inline functions -----------------------------------*/
1638 
1639 
1640 
1641 template <typename number>
1642 inline
1644 {
1645  Assert (cols != 0, ExcNotInitialized());
1646  return cols->rows;
1647 }
1648 
1649 
1650 template <typename number>
1651 inline
1653 {
1654  Assert (cols != 0, ExcNotInitialized());
1655  return cols->cols;
1656 }
1657 
1658 
1659 // Inline the set() and add() functions, since they will be called frequently.
1660 template <typename number>
1661 inline
1662 void
1664  const size_type j,
1665  const number value)
1666 {
1668 
1669  const size_type index = cols->operator()(i, j);
1670 
1671  // it is allowed to set elements of the matrix that are not part of the
1672  // sparsity pattern, if the value to which we set it is zero
1674  {
1676  (value == number()),
1677  ExcInvalidIndex(i, j));
1678  return;
1679  }
1680 
1681  val[index] = value;
1682 }
1683 
1684 
1685 
1686 template <typename number>
1687 template <typename number2>
1688 inline
1689 void
1690 SparseMatrix<number>::set (const std::vector<size_type> &indices,
1691  const FullMatrix<number2> &values,
1692  const bool elide_zero_values)
1693 {
1694  Assert (indices.size() == values.m(),
1695  ExcDimensionMismatch(indices.size(), values.m()));
1696  Assert (values.m() == values.n(), ExcNotQuadratic());
1697 
1698  for (size_type i=0; i<indices.size(); ++i)
1699  set (indices[i], indices.size(), &indices[0], &values(i,0),
1700  elide_zero_values);
1701 }
1702 
1703 
1704 
1705 template <typename number>
1706 template <typename number2>
1707 inline
1708 void
1709 SparseMatrix<number>::set (const std::vector<size_type> &row_indices,
1710  const std::vector<size_type> &col_indices,
1711  const FullMatrix<number2> &values,
1712  const bool elide_zero_values)
1713 {
1714  Assert (row_indices.size() == values.m(),
1715  ExcDimensionMismatch(row_indices.size(), values.m()));
1716  Assert (col_indices.size() == values.n(),
1717  ExcDimensionMismatch(col_indices.size(), values.n()));
1718 
1719  for (size_type i=0; i<row_indices.size(); ++i)
1720  set (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1721  elide_zero_values);
1722 }
1723 
1724 
1725 
1726 template <typename number>
1727 template <typename number2>
1728 inline
1729 void
1731  const std::vector<size_type> &col_indices,
1732  const std::vector<number2> &values,
1733  const bool elide_zero_values)
1734 {
1735  Assert (col_indices.size() == values.size(),
1736  ExcDimensionMismatch(col_indices.size(), values.size()));
1737 
1738  set (row, col_indices.size(), &col_indices[0], &values[0],
1739  elide_zero_values);
1740 }
1741 
1742 
1743 
1744 template <typename number>
1745 inline
1746 void
1748  const size_type j,
1749  const number value)
1750 {
1752 
1753  if (value == number())
1754  return;
1755 
1756  const size_type index = cols->operator()(i, j);
1757 
1758  // it is allowed to add elements to the matrix that are not part of the
1759  // sparsity pattern, if the value to which we set it is zero
1761  {
1763  (value == number()),
1764  ExcInvalidIndex(i, j));
1765  return;
1766  }
1767 
1768  val[index] += value;
1769 }
1770 
1771 
1772 
1773 template <typename number>
1774 template <typename number2>
1775 inline
1776 void
1777 SparseMatrix<number>::add (const std::vector<size_type> &indices,
1778  const FullMatrix<number2> &values,
1779  const bool elide_zero_values)
1780 {
1781  Assert (indices.size() == values.m(),
1782  ExcDimensionMismatch(indices.size(), values.m()));
1783  Assert (values.m() == values.n(), ExcNotQuadratic());
1784 
1785  for (size_type i=0; i<indices.size(); ++i)
1786  add (indices[i], indices.size(), &indices[0], &values(i,0),
1787  elide_zero_values);
1788 }
1789 
1790 
1791 
1792 template <typename number>
1793 template <typename number2>
1794 inline
1795 void
1796 SparseMatrix<number>::add (const std::vector<size_type> &row_indices,
1797  const std::vector<size_type> &col_indices,
1798  const FullMatrix<number2> &values,
1799  const bool elide_zero_values)
1800 {
1801  Assert (row_indices.size() == values.m(),
1802  ExcDimensionMismatch(row_indices.size(), values.m()));
1803  Assert (col_indices.size() == values.n(),
1804  ExcDimensionMismatch(col_indices.size(), values.n()));
1805 
1806  for (size_type i=0; i<row_indices.size(); ++i)
1807  add (row_indices[i], col_indices.size(), &col_indices[0], &values(i,0),
1808  elide_zero_values);
1809 }
1810 
1811 
1812 
1813 template <typename number>
1814 template <typename number2>
1815 inline
1816 void
1818  const std::vector<size_type> &col_indices,
1819  const std::vector<number2> &values,
1820  const bool elide_zero_values)
1821 {
1822  Assert (col_indices.size() == values.size(),
1823  ExcDimensionMismatch(col_indices.size(), values.size()));
1824 
1825  add (row, col_indices.size(), &col_indices[0], &values[0],
1826  elide_zero_values);
1827 }
1828 
1829 
1830 
1831 template <typename number>
1832 inline
1834 SparseMatrix<number>::operator *= (const number factor)
1835 {
1836  Assert (cols != 0, ExcNotInitialized());
1837  Assert (val != 0, ExcNotInitialized());
1838 
1839  number *val_ptr = &val[0];
1840  const number *const end_ptr = &val[cols->n_nonzero_elements()];
1841 
1842  while (val_ptr != end_ptr)
1843  *val_ptr++ *= factor;
1844 
1845  return *this;
1846 }
1847 
1848 
1849 
1850 template <typename number>
1851 inline
1853 SparseMatrix<number>::operator /= (const number factor)
1854 {
1855  Assert (cols != 0, ExcNotInitialized());
1856  Assert (val != 0, ExcNotInitialized());
1857  Assert (factor != number(), ExcDivideByZero());
1858 
1859  const number factor_inv = number(1.) / factor;
1860 
1861  number *val_ptr = &val[0];
1862  const number *const end_ptr = &val[cols->n_nonzero_elements()];
1863 
1864  while (val_ptr != end_ptr)
1865  *val_ptr++ *= factor_inv;
1866 
1867  return *this;
1868 }
1869 
1870 
1871 
1872 template <typename number>
1873 inline
1875  const size_type j) const
1876 {
1877  Assert (cols != 0, ExcNotInitialized());
1878  Assert (cols->operator()(i,j) != SparsityPattern::invalid_entry,
1879  ExcInvalidIndex(i,j));
1880  return val[cols->operator()(i,j)];
1881 }
1882 
1883 
1884 
1885 template <typename number>
1886 inline
1887 number SparseMatrix<number>::el (const size_type i,
1888  const size_type j) const
1889 {
1890  Assert (cols != 0, ExcNotInitialized());
1891  const size_type index = cols->operator()(i,j);
1892 
1894  return val[index];
1895  else
1896  return 0;
1897 }
1898 
1899 
1900 
1901 template <typename number>
1902 inline
1903 number SparseMatrix<number>::diag_element (const size_type i) const
1904 {
1905  Assert (cols != 0, ExcNotInitialized());
1906  Assert (m() == n(), ExcNotQuadratic());
1907  AssertIndexRange(i, m());
1908 
1909  // Use that the first element in each row of a quadratic matrix is the main
1910  // diagonal
1911  return val[cols->rowstart[i]];
1912 }
1913 
1914 
1915 
1916 template <typename number>
1917 inline
1919 {
1920  Assert (cols != 0, ExcNotInitialized());
1921  Assert (m() == n(), ExcNotQuadratic());
1922  AssertIndexRange(i, m());
1923 
1924  // Use that the first element in each row of a quadratic matrix is the main
1925  // diagonal
1926  return val[cols->rowstart[i]];
1927 }
1928 
1929 
1930 
1931 template <typename number>
1932 template <typename ForwardIterator>
1933 void
1934 SparseMatrix<number>::copy_from (const ForwardIterator begin,
1935  const ForwardIterator end)
1936 {
1937  Assert (static_cast<size_type>(std::distance (begin, end)) == m(),
1938  ExcIteratorRange (std::distance (begin, end), m()));
1939 
1940  // for use in the inner loop, we define a typedef to the type of the inner
1941  // iterators
1942  typedef typename std::iterator_traits<ForwardIterator>::value_type::const_iterator inner_iterator;
1943  size_type row=0;
1944  for (ForwardIterator i=begin; i!=end; ++i, ++row)
1945  {
1946  const inner_iterator end_of_row = i->end();
1947  for (inner_iterator j=i->begin(); j!=end_of_row; ++j)
1948  // write entries
1949  set (row, j->first, j->second);
1950  };
1951 }
1952 
1953 
1954 //---------------------------------------------------------------------------
1955 
1956 
1957 namespace SparseMatrixIterators
1958 {
1959  template <typename number>
1960  inline
1962  Accessor (const MatrixType *matrix,
1963  const std::size_t index_within_matrix)
1964  :
1965  SparsityPatternIterators::Accessor (&matrix->get_sparsity_pattern(),
1966  index_within_matrix),
1967  matrix (matrix)
1968  {}
1969 
1970 
1971 
1972  template <typename number>
1973  inline
1974  Accessor<number,true>::
1975  Accessor (const MatrixType *matrix)
1976  :
1977  SparsityPatternIterators::Accessor (&matrix->get_sparsity_pattern()),
1978  matrix (matrix)
1979  {}
1980 
1981 
1982 
1983  template <typename number>
1984  inline
1985  Accessor<number,true>::
1987  :
1988  SparsityPatternIterators::Accessor (a),
1989  matrix (&a.get_matrix())
1990  {}
1991 
1992 
1993 
1994  template <typename number>
1995  inline
1996  number
1997  Accessor<number, true>::value () const
1998  {
1999  AssertIndexRange(index_within_sparsity, matrix->n_nonzero_elements());
2000  return matrix->val[index_within_sparsity];
2001  }
2002 
2003 
2004 
2005  template <typename number>
2006  inline
2007  const typename Accessor<number, true>::MatrixType &
2008  Accessor<number, true>::get_matrix () const
2009  {
2010  return *matrix;
2011  }
2012 
2013 
2014 
2015  template <typename number>
2016  inline
2017  Accessor<number, false>::Reference::Reference (
2018  const Accessor *accessor,
2019  const bool)
2020  :
2021  accessor (accessor)
2022  {}
2023 
2024 
2025  template <typename number>
2026  inline
2027  Accessor<number, false>::Reference::operator number() const
2028  {
2029  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2030  return accessor->matrix->val[accessor->index_within_sparsity];
2031  }
2032 
2033 
2034 
2035  template <typename number>
2036  inline
2037  const typename Accessor<number, false>::Reference &
2038  Accessor<number, false>::Reference::operator = (const number n) const
2039  {
2040  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2041  accessor->matrix->val[accessor->index_within_sparsity] = n;
2042  return *this;
2043  }
2044 
2045 
2046 
2047  template <typename number>
2048  inline
2049  const typename Accessor<number, false>::Reference &
2050  Accessor<number, false>::Reference::operator += (const number n) const
2051  {
2052  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2053  accessor->matrix->val[accessor->index_within_sparsity] += n;
2054  return *this;
2055  }
2056 
2057 
2058 
2059  template <typename number>
2060  inline
2061  const typename Accessor<number, false>::Reference &
2062  Accessor<number, false>::Reference::operator -= (const number n) const
2063  {
2064  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2065  accessor->matrix->val[accessor->index_within_sparsity] -= n;
2066  return *this;
2067  }
2068 
2069 
2070 
2071  template <typename number>
2072  inline
2073  const typename Accessor<number, false>::Reference &
2074  Accessor<number, false>::Reference::operator *= (const number n) const
2075  {
2076  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2077  accessor->matrix->val[accessor->index_within_sparsity] *= n;
2078  return *this;
2079  }
2080 
2081 
2082 
2083  template <typename number>
2084  inline
2085  const typename Accessor<number, false>::Reference &
2086  Accessor<number, false>::Reference::operator /= (const number n) const
2087  {
2088  AssertIndexRange(accessor->index_within_sparsity, accessor->matrix->n_nonzero_elements());
2089  accessor->matrix->val[accessor->index_within_sparsity] /= n;
2090  return *this;
2091  }
2092 
2093 
2094 
2095  template <typename number>
2096  inline
2097  Accessor<number,false>::
2098  Accessor (MatrixType *matrix,
2099  const std::size_t index)
2100  :
2101  SparsityPatternIterators::Accessor (&matrix->get_sparsity_pattern(),
2102  index),
2103  matrix (matrix)
2104  {}
2105 
2106 
2107 
2108  template <typename number>
2109  inline
2110  Accessor<number,false>::
2111  Accessor (MatrixType *matrix)
2112  :
2113  SparsityPatternIterators::Accessor (&matrix->get_sparsity_pattern()),
2114  matrix (matrix)
2115  {}
2116 
2117 
2118 
2119  template <typename number>
2120  inline
2121  typename Accessor<number, false>::Reference
2122  Accessor<number, false>::value() const
2123  {
2124  return Reference(this,true);
2125  }
2126 
2127 
2128 
2129 
2130  template <typename number>
2131  inline
2132  typename Accessor<number, false>::MatrixType &
2133  Accessor<number, false>::get_matrix () const
2134  {
2135  return *matrix;
2136  }
2137 
2138 
2139 
2140  template <typename number, bool Constness>
2141  inline
2142  Iterator<number, Constness>::
2143  Iterator (MatrixType *matrix,
2144  const std::size_t index)
2145  :
2146  accessor(matrix, index)
2147  {}
2148 
2149 
2150 
2151  template <typename number, bool Constness>
2152  inline
2153  Iterator<number, Constness>::
2154  Iterator (MatrixType *matrix)
2155  :
2156  accessor(matrix)
2157  {}
2158 
2159 
2160 
2161  template <typename number, bool Constness>
2162  inline
2163  Iterator<number, Constness>::
2165  :
2166  accessor(*i)
2167  {}
2168 
2169 
2170 
2171  template <typename number, bool Constness>
2172  inline
2173  Iterator<number, Constness> &
2174  Iterator<number,Constness>::operator++ ()
2175  {
2176  accessor.advance ();
2177  return *this;
2178  }
2179 
2180 
2181  template <typename number, bool Constness>
2182  inline
2183  Iterator<number,Constness>
2184  Iterator<number,Constness>::operator++ (int)
2185  {
2186  const Iterator iter = *this;
2187  accessor.advance ();
2188  return iter;
2189  }
2190 
2191 
2192  template <typename number, bool Constness>
2193  inline
2194  const Accessor<number,Constness> &
2195  Iterator<number,Constness>::operator* () const
2196  {
2197  return accessor;
2198  }
2199 
2200 
2201  template <typename number, bool Constness>
2202  inline
2203  const Accessor<number,Constness> *
2204  Iterator<number,Constness>::operator-> () const
2205  {
2206  return &accessor;
2207  }
2208 
2209 
2210  template <typename number, bool Constness>
2211  inline
2212  bool
2213  Iterator<number,Constness>::
2214  operator == (const Iterator &other) const
2215  {
2216  return (accessor == other.accessor);
2217  }
2218 
2219 
2220  template <typename number, bool Constness>
2221  inline
2222  bool
2223  Iterator<number,Constness>::
2224  operator != (const Iterator &other) const
2225  {
2226  return ! (*this == other);
2227  }
2228 
2229 
2230  template <typename number, bool Constness>
2231  inline
2232  bool
2233  Iterator<number,Constness>::
2234  operator < (const Iterator &other) const
2235  {
2236  Assert (&accessor.get_matrix() == &other.accessor.get_matrix(),
2237  ExcInternalError());
2238 
2239  return (accessor < other.accessor);
2240  }
2241 
2242 
2243  template <typename number, bool Constness>
2244  inline
2245  bool
2246  Iterator<number,Constness>::
2247  operator > (const Iterator &other) const
2248  {
2249  return (other < *this);
2250  }
2251 
2252 
2253  template <typename number, bool Constness>
2254  inline
2255  int
2256  Iterator<number,Constness>::
2257  operator - (const Iterator &other) const
2258  {
2259  Assert (&accessor.get_matrix() == &other.accessor.get_matrix(),
2260  ExcInternalError());
2261 
2262  return (*this)->index_within_sparsity - other->index_within_sparsity;
2263  }
2264 
2265 
2266 
2267  template <typename number, bool Constness>
2268  inline
2269  Iterator<number,Constness>
2270  Iterator<number,Constness>::
2271  operator + (const size_type n) const
2272  {
2273  Iterator x = *this;
2274  for (size_type i=0; i<n; ++i)
2275  ++x;
2276 
2277  return x;
2278  }
2279 
2280 }
2281 
2282 
2283 
2284 template <typename number>
2285 inline
2288 {
2289  return const_iterator(this, 0);
2290 }
2291 
2292 
2293 template <typename number>
2294 inline
2297 {
2298  return const_iterator(this);
2299 }
2300 
2301 
2302 template <typename number>
2303 inline
2306 {
2307  return iterator (this, 0);
2308 }
2309 
2310 
2311 template <typename number>
2312 inline
2315 {
2316  return iterator(this, cols->rowstart[cols->rows]);
2317 }
2318 
2319 
2320 template <typename number>
2321 inline
2323 SparseMatrix<number>::begin (const size_type r) const
2324 {
2325  Assert (r<m(), ExcIndexRange(r,0,m()));
2326 
2327  return const_iterator(this, cols->rowstart[r]);
2328 }
2329 
2330 
2331 
2332 template <typename number>
2333 inline
2335 SparseMatrix<number>::end (const size_type r) const
2336 {
2337  Assert (r<m(), ExcIndexRange(r,0,m()));
2338 
2339  return const_iterator(this, cols->rowstart[r+1]);
2340 }
2341 
2342 
2343 
2344 template <typename number>
2345 inline
2347 SparseMatrix<number>::begin (const size_type r)
2348 {
2349  Assert (r<m(), ExcIndexRange(r,0,m()));
2350 
2351  return iterator(this, cols->rowstart[r]);
2352 }
2353 
2354 
2355 
2356 template <typename number>
2357 inline
2359 SparseMatrix<number>::end (const size_type r)
2360 {
2361  Assert (r<m(), ExcIndexRange(r,0,m()));
2362 
2363  return iterator(this, cols->rowstart[r+1]);
2364 }
2365 
2366 
2367 
2368 template <typename number>
2369 template <class StreamType>
2370 inline
2371 void SparseMatrix<number>::print (StreamType &out,
2372  const bool across,
2373  const bool diagonal_first) const
2374 {
2375  Assert (cols != 0, ExcNotInitialized());
2376  Assert (val != 0, ExcNotInitialized());
2377 
2378  bool hanging_diagonal = false;
2379  number diagonal = number();
2380 
2381  for (size_type i=0; i<cols->rows; ++i)
2382  {
2383  for (size_type j=cols->rowstart[i]; j<cols->rowstart[i+1]; ++j)
2384  {
2385  if (!diagonal_first && i == cols->colnums[j])
2386  {
2387  diagonal = val[j];
2388  hanging_diagonal = true;
2389  }
2390  else
2391  {
2392  if (hanging_diagonal && cols->colnums[j]>i)
2393  {
2394  if (across)
2395  out << ' ' << i << ',' << i << ':' << diagonal;
2396  else
2397  out << '(' << i << ',' << i << ") " << diagonal << std::endl;
2398  hanging_diagonal = false;
2399  }
2400  if (across)
2401  out << ' ' << i << ',' << cols->colnums[j] << ':' << val[j];
2402  else
2403  out << "(" << i << "," << cols->colnums[j] << ") " << val[j] << std::endl;
2404  }
2405  }
2406  if (hanging_diagonal)
2407  {
2408  if (across)
2409  out << ' ' << i << ',' << i << ':' << diagonal;
2410  else
2411  out << '(' << i << ',' << i << ") " << diagonal << std::endl;
2412  hanging_diagonal = false;
2413  }
2414  }
2415  if (across)
2416  out << std::endl;
2417 }
2418 
2419 
2420 template <typename number>
2421 inline
2422 void
2424 prepare_add()
2425 {
2426  //nothing to do here
2427 }
2428 
2429 
2430 
2431 template <typename number>
2432 inline
2433 void
2435 prepare_set()
2436 {
2437  //nothing to do here
2438 }
2439 
2440 #endif // DOXYGEN
2441 
2442 
2443 /*---------------------------- sparse_matrix.h ---------------------------*/
2444 
2445 DEAL_II_NAMESPACE_CLOSE
2446 
2447 #endif
2448 /*---------------------------- sparse_matrix.h ---------------------------*/
void PSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const
void Tvmult_add(OutVector &dst, const InVector &src) const
size_type m() const
SparseMatrix & operator/=(const number factor)
void precondition_TSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
real_type l1_norm() const
size_type get_row_length(const size_type row) const
const Accessor< number, Constness > & operator*() const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
static ::ExceptionBase & ExcInvalidIndex(int arg1, int arg2)
numbers::NumberTraits< number >::real_type real_type
int operator-(const Iterator &p) const
void prepare_add()
const_iterator end() const
Accessor< number, Constness >::MatrixType MatrixType
SparseMatrix< number > & operator=(const SparseMatrix< number > &)
std::size_t n_actually_nonzero_elements(const double threshold=0.) const
static const size_type invalid_entry
const SparseMatrix< number > & get_matrix() const
Iterator(MatrixType *matrix, const std::size_t index_within_matrix)
void set(const size_type i, const size_type j, const number value)
void TSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
void block_read(std::istream &in)
const Accessor< number, Constness > * operator->() const
somenumber matrix_norm_square(const Vector< somenumber > &v) const
bool operator!=(const Iterator &) const
SparseMatrixIterators::Iterator< number, false > iterator
void TPSOR(Vector< somenumber > &v, const std::vector< size_type > &permutation, const std::vector< size_type > &inverse_permutation, const number om=1.) const
static ::ExceptionBase & ExcNotInitialized()
void vmult_add(OutVector &dst, const InVector &src) const
number operator()(const size_type i, const size_type j) const
bool empty() const
static ::ExceptionBase & ExcDifferentSparsityPatterns()
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual void reinit(const SparsityPattern &sparsity)
static ::ExceptionBase & ExcDivideByZero()
void mmult(SparseMatrix< numberC > &C, const SparseMatrix< numberB > &B, const Vector< number > &V=Vector< number >(), const bool rebuild_sparsity_pattern=true) const
void SOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
std::size_t n_nonzero_elements() const
std::size_t memory_consumption() const
size_type n() const
void precondition_SOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
void symmetrize()
size_type n() const
bool operator>(const Iterator &) const
unsigned int global_dof_index
Definition: types.h:88
number diag_element(const size_type i) const
Iterator operator+(const size_type n) const
#define Assert(cond, exc)
Definition: exceptions.h:313
bool operator<(const Iterator &) const
void vmult(OutVector &dst, const InVector &src) const
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void print(StreamType &out, const bool across=false, const bool diagonal_first=true) const
types::global_dof_index size_type
Definition: sparse_matrix.h:56
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
number value_type
somenumber matrix_scalar_product(const Vector< somenumber > &u, const Vector< somenumber > &v) const
#define DeclException0(Exception0)
Definition: exceptions.h:541
somenumber residual(Vector< somenumber > &dst, const Vector< somenumber > &x, const Vector< somenumber > &b) const
void Tmmult(SparseMatrix< numberC > &C, const SparseMatrix< numberB > &B, const Vector< number > &V=Vector< number >(), const bool rebuild_sparsity_pattern=true) const
void precondition_SSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1., const std::vector< std::size_t > &pos_right_of_diagonal=std::vector< std::size_t >()) const
static ::ExceptionBase & ExcIteratorRange(int arg1, int arg2)
void prepare_set()
bool operator==(const Iterator &) const
void SOR(Vector< somenumber > &v, const number om=1.) const
void SSOR(Vector< somenumber > &v, const number omega=1.) const
void add(const size_type i, const size_type j, const number value)
const Accessor< number, Constness > & value_type
Accessor< number, Constness > accessor
void Jacobi_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
SparseMatrixIterators::Iterator< number, true > const_iterator
static ::ExceptionBase & ExcNotQuadratic()
types::global_dof_index size_type
void SSOR_step(Vector< somenumber > &v, const Vector< somenumber > &b, const number om=1.) const
number el(const size_type i, const size_type j) const
SmartPointer< const SparsityPattern, SparseMatrix< number > > cols
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
virtual ~SparseMatrix()
size_type m() const
void TSOR(Vector< somenumber > &v, const number om=1.) const
Accessor(const SparsityPattern *matrix, const std::size_t index_within_sparsity)
static const bool zero_addition_can_be_elided
const SparsityPattern & get_sparsity_pattern() const
const_iterator begin() const
SparseMatrix< number > & copy_from(const SparseMatrix< somenumber > &source)
void print_pattern(std::ostream &out, const double threshold=0.) const
real_type frobenius_norm() const
virtual void clear()
#define AssertIsFinite(number)
Definition: exceptions.h:1197
real_type linfty_norm() const
void Tvmult(OutVector &dst, const InVector &src) const
std::size_t max_len
void compress(::VectorOperation::values)
static ::ExceptionBase & ExcSourceEqualsDestination()
SparseMatrix & operator*=(const number factor)
static ::ExceptionBase & ExcInternalError()
void block_write(std::ostream &out) const