Reference documentation for deal.II version 8.5.1
trilinos_sparse_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 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__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 #ifdef DEAL_II_WITH_CXX11
33 #include <deal.II/lac/vector_memory.h>
34 #include <type_traits>
35 #endif // DEAL_II_WITH_CXX11
36 
37 # include <vector>
38 # include <cmath>
39 # include <memory>
40 
41 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
42 # include <Epetra_FECrsMatrix.h>
43 # include <Epetra_Map.h>
44 # include <Epetra_CrsGraph.h>
45 # include <Epetra_MultiVector.h>
46 # include <Epetra_Operator.h>
47 # include <Epetra_Comm.h>
48 # ifdef DEAL_II_WITH_MPI
49 # include <Epetra_MpiComm.h>
50 # include "mpi.h"
51 # else
52 # include <Epetra_SerialComm.h>
53 # endif
54 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
55 
56 class Epetra_Export;
57 
58 DEAL_II_NAMESPACE_OPEN
59 
60 // forward declarations
61 template <typename MatrixType> class BlockMatrixBase;
62 
63 template <typename number> class SparseMatrix;
64 class SparsityPattern;
66 
67 
68 
69 namespace TrilinosWrappers
70 {
71  // forward declarations
72  class SparseMatrix;
73  class SparsityPattern;
74 
79  {
80  // forward declaration
81  template <bool Constness> class Iterator;
82 
87 
92  std::size_t, std::size_t, std::size_t,
93  << "You tried to access row " << arg1
94  << " of a distributed sparsity pattern, "
95  << " but only rows " << arg2 << " through " << arg3
96  << " are stored locally and can be accessed.");
97 
109  {
110  public:
114  typedef ::types::global_dof_index size_type;
115 
120  const size_type row,
121  const size_type index);
122 
126  size_type row() const;
127 
131  size_type index() const;
132 
136  size_type column() const;
137 
138  protected:
150 
155 
161  void visit_present_row ();
162 
175  std_cxx11::shared_ptr<std::vector<size_type> > colnum_cache;
176 
180  std_cxx11::shared_ptr<std::vector<TrilinosScalar> > value_cache;
181  };
182 
193  template <bool Constess>
194  class Accessor : public AccessorBase
195  {
199  TrilinosScalar value() const;
200 
204  TrilinosScalar &value();
205  };
206 
210  template<>
211  class Accessor<true> : public AccessorBase
212  {
213  public:
218  typedef const SparseMatrix MatrixType;
219 
225  const size_type row,
226  const size_type index);
227 
232  template <bool Other>
233  Accessor (const Accessor<Other> &a);
234 
238  TrilinosScalar value() const;
239 
240  private:
244  template <bool> friend class Iterator;
245  };
246 
250  template<>
251  class Accessor<false> : public AccessorBase
252  {
253  class Reference
254  {
255  public:
259  Reference (const Accessor<false> &accessor);
260 
264  operator TrilinosScalar () const;
265 
269  const Reference &operator = (const TrilinosScalar n) const;
270 
274  const Reference &operator += (const TrilinosScalar n) const;
275 
279  const Reference &operator -= (const TrilinosScalar n) const;
280 
284  const Reference &operator *= (const TrilinosScalar n) const;
285 
289  const Reference &operator /= (const TrilinosScalar n) const;
290 
291  private:
296  Accessor &accessor;
297  };
298 
299  public:
305 
311  const size_type row,
312  const size_type index);
313 
317  Reference value() const;
318 
319  private:
323  template <bool> friend class Iterator;
327  friend class Reference;
328  };
329 
344  template <bool Constness>
345  class Iterator
346  {
347  public:
351  typedef ::types::global_dof_index size_type;
352 
358 
363  Iterator (MatrixType *matrix,
364  const size_type row,
365  const size_type index);
366 
370  template <bool Other>
371  Iterator(const Iterator<Other> &other);
372 
377 
382 
386  const Accessor<Constness> &operator* () const;
387 
391  const Accessor<Constness> *operator-> () const;
392 
397  bool operator == (const Iterator<Constness> &) const;
398 
402  bool operator != (const Iterator<Constness> &) const;
403 
409  bool operator < (const Iterator<Constness> &) const;
410 
414  bool operator > (const Iterator<Constness> &) const;
415 
421  << "Attempt to access element " << arg2
422  << " of row " << arg1
423  << " which doesn't have that many elements.");
424 
425  private:
430 
431  template <bool Other> friend class Iterator;
432  };
433 
434  }
435 
436 
498  class SparseMatrix : public Subscriptor
499  {
500  public:
504  typedef ::types::global_dof_index size_type;
505 
513  struct Traits
514  {
519  static const bool zero_addition_can_be_elided = true;
520  };
521 
526 
531 
535  typedef TrilinosScalar value_type;
536 
544  SparseMatrix ();
545 
553  SparseMatrix (const size_type m,
554  const size_type n,
555  const unsigned int n_max_entries_per_row);
556 
564  SparseMatrix (const size_type m,
565  const size_type n,
566  const std::vector<unsigned int> &n_entries_per_row);
567 
571  SparseMatrix (const SparsityPattern &InputSparsityPattern);
572 
576  virtual ~SparseMatrix ();
577 
593  template<typename SparsityPatternType>
594  void reinit (const SparsityPatternType &sparsity_pattern);
595 
608  void reinit (const SparsityPattern &sparsity_pattern);
609 
618  void reinit (const SparseMatrix &sparse_matrix);
619 
640  template <typename number>
641  void reinit (const ::SparseMatrix<number> &dealii_sparse_matrix,
642  const double drop_tolerance=1e-13,
643  const bool copy_values=true,
644  const ::SparsityPattern *use_this_sparsity=0);
645 
651  void reinit (const Epetra_CrsMatrix &input_matrix,
652  const bool copy_values = true);
654 
672  SparseMatrix (const Epetra_Map &parallel_partitioning,
673  const size_type n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
674 
684  SparseMatrix (const Epetra_Map &parallel_partitioning,
685  const std::vector<unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
686 
705  SparseMatrix (const Epetra_Map &row_parallel_partitioning,
706  const Epetra_Map &col_parallel_partitioning,
707  const size_type n_max_entries_per_row = 0) DEAL_II_DEPRECATED;
708 
725  SparseMatrix (const Epetra_Map &row_parallel_partitioning,
726  const Epetra_Map &col_parallel_partitioning,
727  const std::vector<unsigned int> &n_entries_per_row) DEAL_II_DEPRECATED;
728 
755  template<typename SparsityPatternType>
756  void reinit (const Epetra_Map &parallel_partitioning,
757  const SparsityPatternType &sparsity_pattern,
758  const bool exchange_data = false) DEAL_II_DEPRECATED;
759 
774  template<typename SparsityPatternType>
775  void reinit (const Epetra_Map &row_parallel_partitioning,
776  const Epetra_Map &col_parallel_partitioning,
777  const SparsityPatternType &sparsity_pattern,
778  const bool exchange_data = false) DEAL_II_DEPRECATED;
779 
798  template <typename number>
799  void reinit (const Epetra_Map &parallel_partitioning,
800  const ::SparseMatrix<number> &dealii_sparse_matrix,
801  const double drop_tolerance=1e-13,
802  const bool copy_values=true,
803  const ::SparsityPattern *use_this_sparsity=0) DEAL_II_DEPRECATED;
804 
820  template <typename number>
821  void reinit (const Epetra_Map &row_parallel_partitioning,
822  const Epetra_Map &col_parallel_partitioning,
823  const ::SparseMatrix<number> &dealii_sparse_matrix,
824  const double drop_tolerance=1e-13,
825  const bool copy_values=true,
826  const ::SparsityPattern *use_this_sparsity=0) DEAL_II_DEPRECATED;
828 
844  SparseMatrix (const IndexSet &parallel_partitioning,
845  const MPI_Comm &communicator = MPI_COMM_WORLD,
846  const unsigned int n_max_entries_per_row = 0);
847 
855  SparseMatrix (const IndexSet &parallel_partitioning,
856  const MPI_Comm &communicator,
857  const std::vector<unsigned int> &n_entries_per_row);
858 
873  SparseMatrix (const IndexSet &row_parallel_partitioning,
874  const IndexSet &col_parallel_partitioning,
875  const MPI_Comm &communicator = MPI_COMM_WORLD,
876  const size_type n_max_entries_per_row = 0);
877 
892  SparseMatrix (const IndexSet &row_parallel_partitioning,
893  const IndexSet &col_parallel_partitioning,
894  const MPI_Comm &communicator,
895  const std::vector<unsigned int> &n_entries_per_row);
896 
917  template<typename SparsityPatternType>
918  void reinit (const IndexSet &parallel_partitioning,
919  const SparsityPatternType &sparsity_pattern,
920  const MPI_Comm &communicator = MPI_COMM_WORLD,
921  const bool exchange_data = false);
922 
935  template<typename SparsityPatternType>
936  void reinit (const IndexSet &row_parallel_partitioning,
937  const IndexSet &col_parallel_partitioning,
938  const SparsityPatternType &sparsity_pattern,
939  const MPI_Comm &communicator = MPI_COMM_WORLD,
940  const bool exchange_data = false);
941 
958  template <typename number>
959  void reinit (const IndexSet &parallel_partitioning,
960  const ::SparseMatrix<number> &dealii_sparse_matrix,
961  const MPI_Comm &communicator = MPI_COMM_WORLD,
962  const double drop_tolerance=1e-13,
963  const bool copy_values=true,
964  const ::SparsityPattern *use_this_sparsity=0);
965 
979  template <typename number>
980  void reinit (const IndexSet &row_parallel_partitioning,
981  const IndexSet &col_parallel_partitioning,
982  const ::SparseMatrix<number> &dealii_sparse_matrix,
983  const MPI_Comm &communicator = MPI_COMM_WORLD,
984  const double drop_tolerance=1e-13,
985  const bool copy_values=true,
986  const ::SparsityPattern *use_this_sparsity=0);
988 
992 
996  size_type m () const;
997 
1001  size_type n () const;
1002 
1011  unsigned int local_size () const;
1012 
1021  std::pair<size_type, size_type>
1022  local_range () const;
1023 
1028  bool in_local_range (const size_type index) const;
1029 
1034  size_type n_nonzero_elements () const;
1035 
1039  unsigned int row_length (const size_type row) const;
1040 
1047  bool is_compressed () const;
1048 
1054  size_type memory_consumption () const;
1055 
1059  MPI_Comm get_mpi_communicator () const;
1060 
1062 
1066 
1076  SparseMatrix &
1077  operator = (const double d);
1078 
1086  void clear ();
1087 
1115  void compress (::VectorOperation::values operation);
1116 
1138  void set (const size_type i,
1139  const size_type j,
1140  const TrilinosScalar value);
1141 
1174  void set (const std::vector<size_type> &indices,
1175  const FullMatrix<TrilinosScalar> &full_matrix,
1176  const bool elide_zero_values = false);
1177 
1183  void set (const std::vector<size_type> &row_indices,
1184  const std::vector<size_type> &col_indices,
1185  const FullMatrix<TrilinosScalar> &full_matrix,
1186  const bool elide_zero_values = false);
1187 
1215  void set (const size_type row,
1216  const std::vector<size_type> &col_indices,
1217  const std::vector<TrilinosScalar> &values,
1218  const bool elide_zero_values = false);
1219 
1247  void set (const size_type row,
1248  const size_type n_cols,
1249  const size_type *col_indices,
1250  const TrilinosScalar *values,
1251  const bool elide_zero_values = false);
1252 
1262  void add (const size_type i,
1263  const size_type j,
1264  const TrilinosScalar value);
1265 
1284  void add (const std::vector<size_type> &indices,
1285  const FullMatrix<TrilinosScalar> &full_matrix,
1286  const bool elide_zero_values = true);
1287 
1293  void add (const std::vector<size_type> &row_indices,
1294  const std::vector<size_type> &col_indices,
1295  const FullMatrix<TrilinosScalar> &full_matrix,
1296  const bool elide_zero_values = true);
1297 
1311  void add (const size_type row,
1312  const std::vector<size_type> &col_indices,
1313  const std::vector<TrilinosScalar> &values,
1314  const bool elide_zero_values = true);
1315 
1329  void add (const size_type row,
1330  const size_type n_cols,
1331  const size_type *col_indices,
1332  const TrilinosScalar *values,
1333  const bool elide_zero_values = true,
1334  const bool col_indices_are_sorted = false);
1335 
1339  SparseMatrix &operator *= (const TrilinosScalar factor);
1340 
1344  SparseMatrix &operator /= (const TrilinosScalar factor);
1345 
1349  void copy_from (const SparseMatrix &source);
1350 
1358  void add (const TrilinosScalar factor,
1359  const SparseMatrix &matrix);
1360 
1387  void clear_row (const size_type row,
1388  const TrilinosScalar new_diag_value = 0);
1389 
1410  void clear_rows (const std::vector<size_type> &rows,
1411  const TrilinosScalar new_diag_value = 0);
1412 
1422  void transpose ();
1423 
1425 
1429 
1438  TrilinosScalar operator () (const size_type i,
1439  const size_type j) const;
1440 
1457  TrilinosScalar el (const size_type i,
1458  const size_type j) const;
1459 
1466  TrilinosScalar diag_element (const size_type i) const;
1467 
1469 
1473 
1495  template<typename VectorType>
1496  void vmult (VectorType &dst,
1497  const VectorType &src) const;
1498 
1521  template <typename VectorType>
1522  void Tvmult (VectorType &dst,
1523  const VectorType &src) const;
1524 
1548  template<typename VectorType>
1549  void vmult_add (VectorType &dst,
1550  const VectorType &src) const;
1551 
1575  template <typename VectorType>
1576  void Tvmult_add (VectorType &dst,
1577  const VectorType &src) const;
1578 
1603  TrilinosScalar matrix_norm_square (const VectorBase &v) const;
1604 
1623  TrilinosScalar matrix_scalar_product (const VectorBase &u,
1624  const VectorBase &v) const;
1625 
1643  TrilinosScalar residual (VectorBase &dst,
1644  const VectorBase &x,
1645  const VectorBase &b) const;
1646 
1661  void mmult (SparseMatrix &C,
1662  const SparseMatrix &B,
1663  const VectorBase &V = VectorBase()) const;
1664 
1665 
1682  void Tmmult (SparseMatrix &C,
1683  const SparseMatrix &B,
1684  const VectorBase &V = VectorBase()) const;
1685 
1687 
1691 
1699  TrilinosScalar l1_norm () const;
1700 
1709  TrilinosScalar linfty_norm () const;
1710 
1715  TrilinosScalar frobenius_norm () const;
1716 
1718 
1722 
1727  const Epetra_CrsMatrix &trilinos_matrix () const;
1728 
1733  const Epetra_CrsGraph &trilinos_sparsity_pattern () const;
1734 
1742  const Epetra_Map &domain_partitioner () const DEAL_II_DEPRECATED;
1743 
1752  const Epetra_Map &range_partitioner () const DEAL_II_DEPRECATED;
1753 
1761  const Epetra_Map &row_partitioner () const DEAL_II_DEPRECATED;
1762 
1772  const Epetra_Map &col_partitioner () const DEAL_II_DEPRECATED;
1774 
1779 
1785 
1792 
1794 
1799 
1818  const_iterator begin () const;
1819 
1823  iterator begin ();
1824 
1829  const_iterator end () const;
1830 
1834  iterator end ();
1835 
1864  const_iterator begin (const size_type r) const;
1865 
1869  iterator begin (const size_type r);
1870 
1880  const_iterator end (const size_type r) const;
1881 
1885  iterator end (const size_type r);
1886 
1888 
1892 
1898  void write_ascii ();
1899 
1907  void print (std::ostream &out,
1908  const bool write_extended_trilinos_info = false) const;
1909 
1911 
1920  int,
1921  << "An error with error number " << arg1
1922  << " occurred while calling a Trilinos function");
1923 
1929  << "The entry with index <" << arg1 << ',' << arg2
1930  << "> does not exist.");
1931 
1936 
1941 
1947  << "You tried to access element (" << arg1
1948  << "/" << arg2 << ")"
1949  << " of a distributed matrix, but only rows "
1950  << arg3 << " through " << arg4
1951  << " are stored locally and can be accessed.");
1952 
1958  << "You tried to access element (" << arg1
1959  << "/" << arg2 << ")"
1960  << " of a sparse matrix, but it appears to not"
1961  << " exist in the Trilinos sparsity pattern.");
1963 
1964 
1965 
1966  protected:
1967 
1978  void prepare_add();
1979 
1985  void prepare_set();
1986 
1987 
1988 
1989  private:
1993  SparseMatrix (const SparseMatrix &);
1997  SparseMatrix &operator = (const SparseMatrix &);
1998 
2003  std_cxx11::shared_ptr<Epetra_Map> column_space_map;
2004 
2010  std_cxx11::shared_ptr<Epetra_FECrsMatrix> matrix;
2011 
2017  std_cxx11::shared_ptr<Epetra_CrsMatrix> nonlocal_matrix;
2018 
2022  std_cxx11::shared_ptr<Epetra_Export> nonlocal_matrix_exporter;
2023 
2035  Epetra_CombineMode last_action;
2036 
2042 
2047  };
2048 
2049 
2050 #ifdef DEAL_II_WITH_CXX11
2051 
2052  // forwards declarations
2053  class SolverBase;
2054  class PreconditionBase;
2055 
2056  namespace internal
2057  {
2058  namespace
2059  {
2060  inline
2061  void check_vector_map_equality(const Epetra_CrsMatrix &mtrx,
2062  const Epetra_MultiVector &src,
2063  const Epetra_MultiVector &dst,
2064  const bool transpose)
2065  {
2066  if (transpose == false)
2067  {
2068  Assert (src.Map().SameAs(mtrx.DomainMap()) == true,
2069  ExcMessage ("Column map of matrix does not fit with vector map!"));
2070  Assert (dst.Map().SameAs(mtrx.RangeMap()) == true,
2071  ExcMessage ("Row map of matrix does not fit with vector map!"));
2072  }
2073  else
2074  {
2075  Assert (src.Map().SameAs(mtrx.RangeMap()) == true,
2076  ExcMessage ("Column map of matrix does not fit with vector map!"));
2077  Assert (dst.Map().SameAs(mtrx.DomainMap()) == true,
2078  ExcMessage ("Row map of matrix does not fit with vector map!"));
2079  }
2080  (void)mtrx; // removes -Wunused-variable in optimized mode
2081  (void)src;
2082  (void)dst;
2083  }
2084 
2085  inline
2086  void check_vector_map_equality(const Epetra_Operator &op,
2087  const Epetra_MultiVector &src,
2088  const Epetra_MultiVector &dst,
2089  const bool transpose)
2090  {
2091  if (transpose == false)
2092  {
2093  Assert (src.Map().SameAs(op.OperatorDomainMap()) == true,
2094  ExcMessage ("Column map of operator does not fit with vector map!"));
2095  Assert (dst.Map().SameAs(op.OperatorRangeMap()) == true,
2096  ExcMessage ("Row map of operator does not fit with vector map!"));
2097  }
2098  else
2099  {
2100  Assert (src.Map().SameAs(op.OperatorRangeMap()) == true,
2101  ExcMessage ("Column map of operator does not fit with vector map!"));
2102  Assert (dst.Map().SameAs(op.OperatorDomainMap()) == true,
2103  ExcMessage ("Row map of operator does not fit with vector map!"));
2104  }
2105  (void)op; // removes -Wunused-variable in optimized mode
2106  (void)src;
2107  (void)dst;
2108  }
2109  }
2110 
2111  namespace LinearOperator
2112  {
2113 
2134  : public Epetra_Operator
2135  {
2136  public:
2137 
2141  typedef Epetra_MultiVector VectorType;
2142 
2147 
2152 
2157 
2165  TrilinosPayload ();
2166 
2170  TrilinosPayload (const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2172 
2176  TrilinosPayload (const TrilinosWrappers::SparseMatrix &matrix_exemplar,
2177  const TrilinosWrappers::PreconditionBase &preconditioner);
2178 
2182  TrilinosPayload (const TrilinosWrappers::PreconditionBase &preconditioner_exemplar,
2183  const TrilinosWrappers::PreconditionBase &preconditioner);
2184 
2188  TrilinosPayload (const TrilinosPayload &payload);
2189 
2197  TrilinosPayload (const TrilinosPayload &first_op,
2198  const TrilinosPayload &second_op);
2199 
2203  virtual ~TrilinosPayload();
2204 
2208  TrilinosPayload identity_payload () const;
2209 
2213  TrilinosPayload null_payload () const;
2214 
2218  TrilinosPayload transpose_payload () const;
2219 
2236  template <typename Solver, typename Preconditioner>
2237  typename std::enable_if<
2238  std::is_base_of<TrilinosWrappers::SolverBase,Solver>::value &&
2239  std::is_base_of<TrilinosWrappers::PreconditionBase,Preconditioner>::value,
2240  TrilinosPayload>::type
2241  inverse_payload (Solver &, const Preconditioner &) const;
2242 
2259  template <typename Solver, typename Preconditioner>
2260  typename std::enable_if<
2261  !(std::is_base_of<TrilinosWrappers::SolverBase,Solver>::value &&
2262  std::is_base_of<TrilinosWrappers::PreconditionBase,Preconditioner>::value),
2263  TrilinosPayload>::type
2264  inverse_payload (Solver &, const Preconditioner &) const;
2265 
2267 
2272 
2278  IndexSet
2280 
2286  IndexSet
2287  locally_owned_range_indices () const;
2288 
2292  MPI_Comm
2293  get_mpi_communicator () const;
2294 
2301  void
2302  transpose ();
2303 
2311  std::function<void(VectorType &, const VectorType &)> vmult;
2312 
2320  std::function<void(VectorType &, const VectorType &)> Tvmult;
2321 
2329  std::function<void(VectorType &, const VectorType &)> inv_vmult;
2330 
2338  std::function<void(VectorType &, const VectorType &)> inv_Tvmult;
2339 
2341 
2346 
2353  virtual bool
2354  UseTranspose () const;
2355 
2371  virtual int
2372  SetUseTranspose (bool UseTranspose);
2373 
2385  virtual int
2386  Apply(const VectorType &X,
2387  VectorType &Y) const;
2388 
2406  virtual int
2407  ApplyInverse(const VectorType &Y,
2408  VectorType &X) const;
2410 
2415 
2422  virtual const char *
2423  Label () const;
2424 
2432  virtual const Epetra_Comm &
2433  Comm () const;
2434 
2442  virtual const Epetra_Map &
2443  OperatorDomainMap () const;
2444 
2453  virtual const Epetra_Map &
2454  OperatorRangeMap () const;
2456 
2457  private:
2458 
2464 
2469 #ifdef DEAL_II_WITH_MPI
2470  Epetra_MpiComm communicator;
2471 #else
2472  Epetra_SerialComm communicator;
2473 #endif
2474 
2479  Epetra_Map domain_map;
2480 
2485  Epetra_Map range_map;
2486 
2495  virtual bool
2496  HasNormInf () const;
2497 
2505  virtual double
2506  NormInf () const;
2507  };
2508 
2513  TrilinosPayload operator+(const TrilinosPayload &first_op,
2514  const TrilinosPayload &second_op);
2515 
2520  TrilinosPayload operator*(const TrilinosPayload &first_op,
2521  const TrilinosPayload &second_op);
2522 
2523  } /* namespace LinearOperator */
2524  } /* namespace internal */
2525 
2526 #endif // DEAL_II_WITH_CXX11
2527 
2528 
2529 // -------------------------- inline and template functions ----------------------
2530 
2531 #ifndef DOXYGEN
2532 
2533  namespace SparseMatrixIterators
2534  {
2535  inline
2536  AccessorBase::AccessorBase(SparseMatrix *matrix, size_type row, size_type index)
2537  :
2538  matrix(matrix),
2539  a_row(row),
2540  a_index(index)
2541  {
2542  visit_present_row ();
2543  }
2544 
2545 
2546  inline
2547  AccessorBase::size_type
2548  AccessorBase::row() const
2549  {
2550  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2551  return a_row;
2552  }
2553 
2554 
2555  inline
2556  AccessorBase::size_type
2557  AccessorBase::column() const
2558  {
2559  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2560  return (*colnum_cache)[a_index];
2561  }
2562 
2563 
2564  inline
2565  AccessorBase::size_type
2566  AccessorBase::index() const
2567  {
2568  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2569  return a_index;
2570  }
2571 
2572 
2573  inline
2574  Accessor<true>::Accessor (MatrixType *matrix,
2575  const size_type row,
2576  const size_type index)
2577  :
2578  AccessorBase(const_cast<SparseMatrix *>(matrix), row, index)
2579  {}
2580 
2581 
2582  template <bool Other>
2583  inline
2584  Accessor<true>::Accessor(const Accessor<Other> &other)
2585  :
2586  AccessorBase(other)
2587  {}
2588 
2589 
2590  inline
2591  TrilinosScalar
2592  Accessor<true>::value() const
2593  {
2594  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2595  return (*value_cache)[a_index];
2596  }
2597 
2598 
2599  inline
2600  Accessor<false>::Reference::Reference (
2601  const Accessor<false> &acc)
2602  :
2603  accessor(const_cast<Accessor<false>&>(acc))
2604  {}
2605 
2606 
2607  inline
2608  Accessor<false>::Reference::operator TrilinosScalar () const
2609  {
2610  return (*accessor.value_cache)[accessor.a_index];
2611  }
2612 
2613  inline
2614  const Accessor<false>::Reference &
2615  Accessor<false>::Reference::operator = (const TrilinosScalar n) const
2616  {
2617  (*accessor.value_cache)[accessor.a_index] = n;
2618  accessor.matrix->set(accessor.row(), accessor.column(),
2619  static_cast<TrilinosScalar>(*this));
2620  return *this;
2621  }
2622 
2623 
2624  inline
2625  const Accessor<false>::Reference &
2626  Accessor<false>::Reference::operator += (const TrilinosScalar n) const
2627  {
2628  (*accessor.value_cache)[accessor.a_index] += n;
2629  accessor.matrix->set(accessor.row(), accessor.column(),
2630  static_cast<TrilinosScalar>(*this));
2631  return *this;
2632  }
2633 
2634 
2635  inline
2636  const Accessor<false>::Reference &
2637  Accessor<false>::Reference::operator -= (const TrilinosScalar n) const
2638  {
2639  (*accessor.value_cache)[accessor.a_index] -= n;
2640  accessor.matrix->set(accessor.row(), accessor.column(),
2641  static_cast<TrilinosScalar>(*this));
2642  return *this;
2643  }
2644 
2645 
2646  inline
2647  const Accessor<false>::Reference &
2648  Accessor<false>::Reference::operator *= (const TrilinosScalar n) const
2649  {
2650  (*accessor.value_cache)[accessor.a_index] *= n;
2651  accessor.matrix->set(accessor.row(), accessor.column(),
2652  static_cast<TrilinosScalar>(*this));
2653  return *this;
2654  }
2655 
2656 
2657  inline
2658  const Accessor<false>::Reference &
2659  Accessor<false>::Reference::operator /= (const TrilinosScalar n) const
2660  {
2661  (*accessor.value_cache)[accessor.a_index] /= n;
2662  accessor.matrix->set(accessor.row(), accessor.column(),
2663  static_cast<TrilinosScalar>(*this));
2664  return *this;
2665  }
2666 
2667 
2668  inline
2669  Accessor<false>::Accessor (MatrixType *matrix,
2670  const size_type row,
2671  const size_type index)
2672  :
2673  AccessorBase(matrix, row, index)
2674  {}
2675 
2676 
2677  inline
2678  Accessor<false>::Reference
2679  Accessor<false>::value() const
2680  {
2681  Assert (a_row < matrix->m(), ExcBeyondEndOfMatrix());
2682  return Reference(*this);
2683  }
2684 
2685 
2686 
2687  template <bool Constness>
2688  inline
2689  Iterator<Constness>::Iterator(MatrixType *matrix,
2690  const size_type row,
2691  const size_type index)
2692  :
2693  accessor(matrix, row, index)
2694  {}
2695 
2696 
2697  template <bool Constness>
2698  template <bool Other>
2699  inline
2700  Iterator<Constness>::Iterator(const Iterator<Other> &other)
2701  :
2702  accessor(other.accessor)
2703  {}
2704 
2705 
2706  template <bool Constness>
2707  inline
2708  Iterator<Constness> &
2709  Iterator<Constness>::operator++ ()
2710  {
2711  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
2712 
2713  ++accessor.a_index;
2714 
2715  // If at end of line: do one
2716  // step, then cycle until we
2717  // find a row with a nonzero
2718  // number of entries.
2719  if (accessor.a_index >= accessor.colnum_cache->size())
2720  {
2721  accessor.a_index = 0;
2722  ++accessor.a_row;
2723 
2724  while ((accessor.a_row < accessor.matrix->m())
2725  &&
2726  ((accessor.matrix->in_local_range (accessor.a_row) == false)
2727  ||
2728  (accessor.matrix->row_length(accessor.a_row) == 0)))
2729  ++accessor.a_row;
2730 
2731  accessor.visit_present_row();
2732  }
2733  return *this;
2734  }
2735 
2736 
2737  template <bool Constness>
2738  inline
2739  Iterator<Constness>
2740  Iterator<Constness>::operator++ (int)
2741  {
2742  const Iterator<Constness> old_state = *this;
2743  ++(*this);
2744  return old_state;
2745  }
2746 
2747 
2748 
2749  template <bool Constness>
2750  inline
2751  const Accessor<Constness> &
2752  Iterator<Constness>::operator* () const
2753  {
2754  return accessor;
2755  }
2756 
2757 
2758 
2759  template <bool Constness>
2760  inline
2761  const Accessor<Constness> *
2762  Iterator<Constness>::operator-> () const
2763  {
2764  return &accessor;
2765  }
2766 
2767 
2768 
2769  template <bool Constness>
2770  inline
2771  bool
2772  Iterator<Constness>::operator == (const Iterator<Constness> &other) const
2773  {
2774  return (accessor.a_row == other.accessor.a_row &&
2775  accessor.a_index == other.accessor.a_index);
2776  }
2777 
2778 
2779 
2780  template <bool Constness>
2781  inline
2782  bool
2783  Iterator<Constness>::operator != (const Iterator<Constness> &other) const
2784  {
2785  return ! (*this == other);
2786  }
2787 
2788 
2789 
2790  template <bool Constness>
2791  inline
2792  bool
2793  Iterator<Constness>::operator < (const Iterator<Constness> &other) const
2794  {
2795  return (accessor.row() < other.accessor.row() ||
2796  (accessor.row() == other.accessor.row() &&
2797  accessor.index() < other.accessor.index()));
2798  }
2799 
2800 
2801  template <bool Constness>
2802  inline
2803  bool
2804  Iterator<Constness>::operator > (const Iterator<Constness> &other) const
2805  {
2806  return (other < *this);
2807  }
2808 
2809  }
2810 
2811 
2812 
2813  inline
2815  SparseMatrix::begin() const
2816  {
2817  return begin(0);
2818  }
2819 
2820 
2821 
2822  inline
2824  SparseMatrix::end() const
2825  {
2826  return const_iterator(this, m(), 0);
2827  }
2828 
2829 
2830 
2831  inline
2833  SparseMatrix::begin(const size_type r) const
2834  {
2835  Assert (r < m(), ExcIndexRange(r, 0, m()));
2836  if (in_local_range (r)
2837  &&
2838  (row_length(r) > 0))
2839  return const_iterator(this, r, 0);
2840  else
2841  return end (r);
2842  }
2843 
2844 
2845 
2846  inline
2848  SparseMatrix::end(const size_type r) const
2849  {
2850  Assert (r < m(), ExcIndexRange(r, 0, m()));
2851 
2852  // place the iterator on the first entry
2853  // past this line, or at the end of the
2854  // matrix
2855  for (size_type i=r+1; i<m(); ++i)
2856  if (in_local_range (i)
2857  &&
2858  (row_length(i) > 0))
2859  return const_iterator(this, i, 0);
2860 
2861  // if there is no such line, then take the
2862  // end iterator of the matrix
2863  return end();
2864  }
2865 
2866 
2867 
2868  inline
2871  {
2872  return begin(0);
2873  }
2874 
2875 
2876 
2877  inline
2880  {
2881  return iterator(this, m(), 0);
2882  }
2883 
2884 
2885 
2886  inline
2888  SparseMatrix::begin(const size_type r)
2889  {
2890  Assert (r < m(), ExcIndexRange(r, 0, m()));
2891  if (in_local_range (r)
2892  &&
2893  (row_length(r) > 0))
2894  return iterator(this, r, 0);
2895  else
2896  return end (r);
2897  }
2898 
2899 
2900 
2901  inline
2903  SparseMatrix::end(const size_type r)
2904  {
2905  Assert (r < m(), ExcIndexRange(r, 0, m()));
2906 
2907  // place the iterator on the first entry
2908  // past this line, or at the end of the
2909  // matrix
2910  for (size_type i=r+1; i<m(); ++i)
2911  if (in_local_range (i)
2912  &&
2913  (row_length(i) > 0))
2914  return iterator(this, i, 0);
2915 
2916  // if there is no such line, then take the
2917  // end iterator of the matrix
2918  return end();
2919  }
2920 
2921 
2922 
2923  inline
2924  bool
2925  SparseMatrix::in_local_range (const size_type index) const
2926  {
2927  TrilinosWrappers::types::int_type begin, end;
2928 #ifndef DEAL_II_WITH_64BIT_INDICES
2929  begin = matrix->RowMap().MinMyGID();
2930  end = matrix->RowMap().MaxMyGID()+1;
2931 #else
2932  begin = matrix->RowMap().MinMyGID64();
2933  end = matrix->RowMap().MaxMyGID64()+1;
2934 #endif
2935 
2936  return ((index >= static_cast<size_type>(begin)) &&
2937  (index < static_cast<size_type>(end)));
2938  }
2939 
2940 
2941 
2942  inline
2943  bool
2945  {
2946  return compressed;
2947  }
2948 
2949 
2950 
2951  // Inline the set() and add() functions, since they will be called
2952  // frequently, and the compiler can optimize away some unnecessary loops
2953  // when the sizes are given at compile time.
2954  inline
2955  void
2956  SparseMatrix::set (const size_type i,
2957  const size_type j,
2958  const TrilinosScalar value)
2959  {
2960 
2961  AssertIsFinite(value);
2962 
2963  set (i, 1, &j, &value, false);
2964  }
2965 
2966 
2967 
2968  inline
2969  void
2970  SparseMatrix::set (const std::vector<size_type> &indices,
2971  const FullMatrix<TrilinosScalar> &values,
2972  const bool elide_zero_values)
2973  {
2974  Assert (indices.size() == values.m(),
2975  ExcDimensionMismatch(indices.size(), values.m()));
2976  Assert (values.m() == values.n(), ExcNotQuadratic());
2977 
2978  for (size_type i=0; i<indices.size(); ++i)
2979  set (indices[i], indices.size(), &indices[0], &values(i,0),
2980  elide_zero_values);
2981  }
2982 
2983 
2984 
2985  inline
2986  void
2987  SparseMatrix::add (const size_type i,
2988  const size_type j,
2989  const TrilinosScalar value)
2990  {
2991  AssertIsFinite(value);
2992 
2993  if (value == 0)
2994  {
2995  // we have to check after Insert/Add in any case to be consistent
2996  // with the MPI communication model (see the comments in the
2997  // documentation of TrilinosWrappers::Vector), but we can save some
2998  // work if the addend is zero. However, these actions are done in case
2999  // we pass on to the other function.
3000 
3001  // TODO: fix this (do not run compress here, but fail)
3002  if (last_action == Insert)
3003  {
3004  int ierr;
3005  ierr = matrix->GlobalAssemble(*column_space_map,
3006  matrix->RowMap(), false);
3007 
3008  Assert (ierr == 0, ExcTrilinosError(ierr));
3009  (void)ierr; // removes -Wunused-but-set-variable in optimized mode
3010  }
3011 
3012  last_action = Add;
3013 
3014  return;
3015  }
3016  else
3017  add (i, 1, &j, &value, false);
3018  }
3019 
3020 
3021 
3022  // inline "simple" functions that are called frequently and do only involve
3023  // a call to some Trilinos function.
3024  inline
3026  SparseMatrix::m () const
3027  {
3028 #ifndef DEAL_II_WITH_64BIT_INDICES
3029  return matrix->NumGlobalRows();
3030 #else
3031  return matrix->NumGlobalRows64();
3032 #endif
3033  }
3034 
3035 
3036 
3037  inline
3039  SparseMatrix::n () const
3040  {
3041  // If the matrix structure has not been fixed (i.e., we did not have a
3042  // sparsity pattern), it does not know about the number of columns so we
3043  // must always take this from the additional column space map
3044  Assert(column_space_map.get() != 0, ExcInternalError());
3045 #ifndef DEAL_II_WITH_64BIT_INDICES
3046  return column_space_map->NumGlobalElements();
3047 #else
3048  return column_space_map->NumGlobalElements64();
3049 #endif
3050  }
3051 
3052 
3053 
3054  inline
3055  unsigned int
3056  SparseMatrix::local_size () const
3057  {
3058  return matrix -> NumMyRows();
3059  }
3060 
3061 
3062 
3063  inline
3064  std::pair<SparseMatrix::size_type, SparseMatrix::size_type>
3065  SparseMatrix::local_range () const
3066  {
3067  size_type begin, end;
3068 #ifndef DEAL_II_WITH_64BIT_INDICES
3069  begin = matrix->RowMap().MinMyGID();
3070  end = matrix->RowMap().MaxMyGID()+1;
3071 #else
3072  begin = matrix->RowMap().MinMyGID64();
3073  end = matrix->RowMap().MaxMyGID64()+1;
3074 #endif
3075 
3076  return std::make_pair (begin, end);
3077  }
3078 
3079 
3080 
3081  inline
3084  {
3085 #ifndef DEAL_II_WITH_64BIT_INDICES
3086  return matrix->NumGlobalNonzeros();
3087 #else
3088  return matrix->NumGlobalNonzeros64();
3089 #endif
3090  }
3091 
3092 
3093 
3094  template <typename SparsityPatternType>
3095  inline
3096  void SparseMatrix::reinit (const IndexSet &parallel_partitioning,
3097  const SparsityPatternType &sparsity_pattern,
3098  const MPI_Comm &communicator,
3099  const bool exchange_data)
3100  {
3101  reinit (parallel_partitioning, parallel_partitioning,
3102  sparsity_pattern, communicator, exchange_data);
3103  }
3104 
3105 
3106 
3107  template <typename number>
3108  inline
3109  void SparseMatrix::reinit (const IndexSet &parallel_partitioning,
3110  const ::SparseMatrix<number> &sparse_matrix,
3111  const MPI_Comm &communicator,
3112  const double drop_tolerance,
3113  const bool copy_values,
3114  const ::SparsityPattern *use_this_sparsity)
3115  {
3116  Epetra_Map map = parallel_partitioning.make_trilinos_map (communicator, false);
3117  reinit (parallel_partitioning, parallel_partitioning, sparse_matrix,
3118  drop_tolerance, copy_values, use_this_sparsity);
3119  }
3120 
3121 
3122 
3123  inline
3124  const Epetra_CrsMatrix &
3126  {
3127  return static_cast<const Epetra_CrsMatrix &>(*matrix);
3128  }
3129 
3130 
3131 
3132  inline
3133  const Epetra_CrsGraph &
3135  {
3136  return matrix->Graph();
3137  }
3138 
3139 
3140 
3141  inline
3142  IndexSet
3144  {
3145  return IndexSet(matrix->DomainMap());
3146  }
3147 
3148 
3149 
3150  inline
3151  IndexSet
3153  {
3154  return IndexSet(matrix->RangeMap());
3155  }
3156 
3157 
3158 
3159  inline
3160  void
3162  {
3163  //nothing to do here
3164  }
3165 
3166 
3167 
3168  inline
3169  void
3171  {
3172  //nothing to do here
3173  }
3174 
3175 
3176 #ifdef DEAL_II_WITH_CXX11
3177  namespace internal
3178  {
3179  namespace LinearOperator
3180  {
3181  template <typename Solver, typename Preconditioner>
3182  typename std::enable_if<
3183  std::is_base_of<TrilinosWrappers::SolverBase,Solver>::value &&
3184  std::is_base_of<TrilinosWrappers::PreconditionBase,Preconditioner>::value,
3185  TrilinosPayload>::type
3186  TrilinosPayload::inverse_payload (
3187  Solver &solver,
3188  const Preconditioner &preconditioner) const
3189  {
3190  const auto &payload = *this;
3191 
3192  TrilinosPayload return_op(payload);
3193 
3194  // Capture by copy so the payloads are always valid
3195 
3196  return_op.inv_vmult = [payload,&solver,&preconditioner](
3197  TrilinosPayload::Domain &tril_dst, const TrilinosPayload::Range &tril_src
3198  )
3199  {
3200  // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3201  // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3203  internal::check_vector_map_equality(payload,
3204  tril_src, tril_dst,
3205  !payload.UseTranspose());
3206  solver.solve(payload, tril_dst, tril_src, preconditioner);
3207  };
3208 
3209  return_op.inv_Tvmult = [payload,&solver,&preconditioner](
3210  TrilinosPayload::Range &tril_dst, const TrilinosPayload::Domain &tril_src
3211  )
3212  {
3213  // Duplicated from TrilinosWrappers::PreconditionBase::vmult
3214  // as well as from TrilinosWrappers::SparseMatrix::Tvmult
3216  internal::check_vector_map_equality(payload,
3217  tril_src, tril_dst,
3218  payload.UseTranspose());
3219 
3220  const_cast<TrilinosPayload &>(payload).transpose();
3221  solver.solve(payload, tril_dst, tril_src, preconditioner);
3222  const_cast<TrilinosPayload &>(payload).transpose();
3223  };
3224 
3225  // If the input operator is already setup for transpose operations, then
3226  // we must do similar with its inverse.
3227  if (return_op.UseTranspose() == true)
3228  std::swap(return_op.inv_vmult,
3229  return_op.inv_Tvmult);
3230 
3231  return return_op;
3232  }
3233 
3234  template <typename Solver, typename Preconditioner>
3235  typename std::enable_if<
3236  !(std::is_base_of<TrilinosWrappers::SolverBase,Solver>::value &&
3237  std::is_base_of<TrilinosWrappers::PreconditionBase,Preconditioner>::value),
3238  TrilinosPayload>::type
3239  TrilinosPayload::inverse_payload (
3240  Solver &,
3241  const Preconditioner &) const
3242  {
3243  TrilinosPayload return_op(*this);
3244 
3245  return_op.inv_vmult = [](TrilinosPayload::Domain &,
3246  const TrilinosPayload::Range &)
3247  {
3248  AssertThrow(false,
3249  ExcMessage("Payload inv_vmult disabled because of "
3250  "incompatible solver/preconditioner choice."));
3251  };
3252 
3253  return_op.inv_Tvmult = [](TrilinosPayload::Range &,
3254  const TrilinosPayload::Domain &)
3255  {
3256  AssertThrow(false,
3257  ExcMessage("Payload inv_vmult disabled because of "
3258  "incompatible solver/preconditioner choice."));
3259  };
3260 
3261  return return_op;
3262  }
3263  } // namespace LinearOperator
3264  } // namespace internal
3265 #endif // DEAL_II_WITH_CXX11
3266 
3267 #endif // DOXYGEN
3268 
3269 } /* namespace TrilinosWrappers */
3270 
3271 
3272 DEAL_II_NAMESPACE_CLOSE
3273 
3274 
3275 #endif // DEAL_II_WITH_TRILINOS
3276 
3277 
3278 /*----------------------- trilinos_sparse_matrix.h --------------------*/
3279 
3280 #endif
3281 /*----------------------- trilinos_sparse_matrix.h --------------------*/
std_cxx11::shared_ptr< std::vector< TrilinosScalar > > value_cache
TrilinosScalar diag_element(const size_type i) const
static ::ExceptionBase & ExcTrilinosError(int arg1)
TrilinosScalar el(const size_type i, const size_type j) 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:576
SparseMatrixIterators::Iterator< false > iterator
static ::ExceptionBase & ExcSourceEqualsDestination()
static ::ExceptionBase & ExcBeyondEndOfMatrix()
bool operator==(const Iterator< Constness > &) const
void clear_row(const size_type row, const TrilinosScalar new_diag_value=0)
void vmult_add(VectorType &dst, const VectorType &src) const
::types::global_dof_index size_type
TrilinosScalar residual(VectorBase &dst, const VectorBase &x, const VectorBase &b) const
STL namespace.
static ::ExceptionBase & ExcInvalidIndex(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcAccessToNonlocalRow(std::size_t arg1, std::size_t arg2, std::size_t arg3)
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
static ::ExceptionBase & ExcAccessToNonPresentElement(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcTrilinosError(int arg1)
const Epetra_Map & col_partitioner() const 1
void clear_rows(const std::vector< size_type > &rows, const TrilinosScalar new_diag_value=0)
void Tvmult(VectorType &dst, const VectorType &src) const
std_cxx11::shared_ptr< Epetra_FECrsMatrix > matrix
Epetra_Map make_trilinos_map(const MPI_Comm &communicator=MPI_COMM_WORLD, const bool overlapping=false) const
Definition: index_set.cc:535
bool in_local_range(const size_type index) const
static ::ExceptionBase & ExcMessage(std::string arg1)
const Accessor< Constness > & operator*() const
std_cxx11::shared_ptr< Epetra_CrsMatrix > nonlocal_matrix
size_type n() const
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
void add(const size_type i, const size_type j, const TrilinosScalar value)
std::function< void(VectorType &, const VectorType &)> Tvmult
#define Assert(cond, exc)
Definition: exceptions.h:313
std_cxx11::shared_ptr< Epetra_Map > column_space_map
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void Tmmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
void print(std::ostream &out, const bool write_extended_trilinos_info=false) const
bool operator!=(const Iterator< Constness > &) const
#define DeclException0(Exception0)
Definition: exceptions.h:541
void mmult(SparseMatrix &C, const SparseMatrix &B, const VectorBase &V=VectorBase()) const
const Accessor< Constness > * operator->() const
const Epetra_CrsMatrix & trilinos_matrix() const
AccessorBase(SparseMatrix *matrix, const size_type row, const size_type index)
static ::ExceptionBase & ExcInvalidIndexWithinRow(size_type arg1, size_type arg2)
const Epetra_Map & range_partitioner() const 1
void vmult(VectorType &dst, const VectorType &src) const
IndexSet locally_owned_range_indices() const
std::function< void(VectorType &, const VectorType &)> inv_vmult
size_type n_nonzero_elements() const
static ::ExceptionBase & ExcIteratorPastEnd()
const_iterator begin() const
static ::ExceptionBase & ExcNotQuadratic()
unsigned int row_length(const size_type row) const
void Tvmult_add(VectorType &dst, const VectorType &src) const
const Epetra_Map & domain_partitioner() const 1
SparseMatrixIterators::Iterator< true > const_iterator
Definition: solver.h:325
types::global_dof_index size_type
void reinit(const SparsityPatternType &sparsity_pattern)
TrilinosScalar matrix_norm_square(const VectorBase &v) const
std::function< void(VectorType &, const VectorType &)> vmult
void copy_from(const SparseMatrix &source)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
std::pair< size_type, size_type > local_range() const
TrilinosScalar matrix_scalar_product(const VectorBase &u, const VectorBase &v) const
static ::ExceptionBase & ExcMatrixNotCompressed()
IndexSet locally_owned_domain_indices() const
const Epetra_CrsGraph & trilinos_sparsity_pattern() const
#define DeclException3(Exception3, type1, type2, type3, outsequence)
Definition: exceptions.h:588
void compress(::VectorOperation::values operation)
static ::ExceptionBase & ExcAccessToNonLocalElement(size_type arg1, size_type arg2, size_type arg3, size_type arg4)
void set(const size_type i, const size_type j, const TrilinosScalar value)
#define AssertIsFinite(number)
Definition: exceptions.h:1197
unsigned int local_size() const
std::function< void(VectorType &, const VectorType &)> inv_Tvmult
static ::ExceptionBase & ExcInternalError()
const Epetra_Map & row_partitioner() const 1
std_cxx11::shared_ptr< Epetra_Export > nonlocal_matrix_exporter