Reference documentation for deal.II version 8.5.1
trilinos_precondition.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_precondition_h
17 #define dealii__trilinos_precondition_h
18 
19 
20 #include <deal.II/base/config.h>
21 
22 #ifdef DEAL_II_WITH_TRILINOS
23 
24 # include <deal.II/base/subscriptor.h>
25 # include <deal.II/base/std_cxx11/shared_ptr.h>
26 
27 # include <deal.II/lac/trilinos_vector_base.h>
28 # include <deal.II/lac/la_parallel_vector.h>
29 
30 DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
31 # ifdef DEAL_II_WITH_MPI
32 # include <Epetra_MpiComm.h>
33 # else
34 # include <Epetra_SerialComm.h>
35 # endif
36 # include <Epetra_Map.h>
37 
38 # include <Teuchos_ParameterList.hpp>
39 # include <Epetra_RowMatrix.h>
40 # include <Epetra_Vector.h>
41 DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
42 
43 // forward declarations
44 class Ifpack_Preconditioner;
45 class Ifpack_Chebyshev;
46 namespace ML_Epetra
47 {
48  class MultiLevelPreconditioner;
49 }
50 
51 
52 DEAL_II_NAMESPACE_OPEN
53 
54 // forward declarations
55 template <typename number> class SparseMatrix;
56 template <typename number> class Vector;
57 class SparsityPattern;
58 
63 namespace TrilinosWrappers
64 {
65  // forward declarations
66  class SparseMatrix;
67  class BlockSparseMatrix;
68  class SolverBase;
69 
79  {
80  public:
84  typedef ::types::global_dof_index size_type;
85 
91  {};
92 
99 
104 
109 
114  void clear ();
115 
119  MPI_Comm get_mpi_communicator () const;
120 
130  void transpose ();
131 
135  virtual void vmult (VectorBase &dst,
136  const VectorBase &src) const;
137 
141  virtual void Tvmult (VectorBase &dst,
142  const VectorBase &src) const;
143 
148  virtual void vmult (::Vector<double> &dst,
149  const ::Vector<double> &src) const;
150 
155  virtual void Tvmult (::Vector<double> &dst,
156  const ::Vector<double> &src) const;
157 
163  const ::LinearAlgebra::distributed::Vector<double> &src) const;
164 
170  const ::LinearAlgebra::distributed::Vector<double> &src) const;
171 
181  Epetra_Operator &trilinos_operator () const;
183 
188 
194 
201 
203 
213  std::string,
214  << "The sparse matrix the preconditioner is based on "
215  << "uses a map that is not compatible to the one in vector "
216  << arg1
217  << ". Check preconditioner and matrix setup.");
219 
220  friend class SolverBase;
221 
222  protected:
227  std_cxx11::shared_ptr<Epetra_Operator> preconditioner;
228 
233 #ifdef DEAL_II_WITH_MPI
234  Epetra_MpiComm communicator;
235 #else
236  Epetra_SerialComm communicator;
237 #endif
238 
243  std_cxx11::shared_ptr<Epetra_Map> vector_distributor;
244  };
245 
246 
264  {
265  public:
266 
279  {
284  AdditionalData (const double omega = 1,
285  const double min_diagonal = 0,
286  const unsigned int n_sweeps = 1);
287 
291  double omega;
292 
300  double min_diagonal;
301 
306  unsigned int n_sweeps;
307  };
308 
313  void initialize (const SparseMatrix &matrix,
314  const AdditionalData &additional_data = AdditionalData());
315  };
316 
317 
318 
319 
347  {
348  public:
349 
364  {
371  AdditionalData (const double omega = 1,
372  const double min_diagonal = 0,
373  const unsigned int overlap = 0,
374  const unsigned int n_sweeps = 1);
375 
380  double omega;
381 
389  double min_diagonal;
390 
395  unsigned int overlap;
396 
401  unsigned int n_sweeps;
402  };
403 
409  void initialize (const SparseMatrix &matrix,
410  const AdditionalData &additional_data = AdditionalData());
411  };
412 
413 
414 
415 
443  {
444  public:
445 
460  {
467  AdditionalData (const double omega = 1,
468  const double min_diagonal = 0,
469  const unsigned int overlap = 0,
470  const unsigned int n_sweeps = 1);
471 
476  double omega;
477 
485  double min_diagonal;
486 
491  unsigned int overlap;
492 
497  unsigned int n_sweeps;
498  };
499 
505  void initialize (const SparseMatrix &matrix,
506  const AdditionalData &additional_data = AdditionalData());
507  };
508 
509 
510 
529  {
530  public:
531 
550  {
556  AdditionalData (const unsigned int block_size = 1,
557  const std::string &block_creation_type = "linear",
558  const double omega = 1,
559  const double min_diagonal = 0,
560  const unsigned int n_sweeps = 1);
561 
565  unsigned int block_size;
566 
574  std::string block_creation_type;
575 
580  double omega;
581 
589  double min_diagonal;
590 
595  unsigned int n_sweeps;
596  };
597 
602  void initialize (const SparseMatrix &matrix,
603  const AdditionalData &additional_data = AdditionalData());
604  };
605 
606 
607 
608 
632  {
633  public:
634 
653  {
661  AdditionalData (const unsigned int block_size = 1,
662  const std::string &block_creation_type = "linear",
663  const double omega = 1,
664  const double min_diagonal = 0,
665  const unsigned int overlap = 0,
666  const unsigned int n_sweeps = 1);
667 
671  unsigned int block_size;
672 
680  std::string block_creation_type;
681 
686  double omega;
687 
695  double min_diagonal;
696 
701  unsigned int overlap;
702 
707  unsigned int n_sweeps;
708  };
709 
715  void initialize (const SparseMatrix &matrix,
716  const AdditionalData &additional_data = AdditionalData());
717  };
718 
719 
720 
721 
745  {
746  public:
747 
766  {
774  AdditionalData (const unsigned int block_size = 1,
775  const std::string &block_creation_type = "linear",
776  const double omega = 1,
777  const double min_diagonal = 0,
778  const unsigned int overlap = 0,
779  const unsigned int n_sweeps = 1);
780 
784  unsigned int block_size;
785 
793  std::string block_creation_type;
794 
799  double omega;
800 
808  double min_diagonal;
809 
814  unsigned int overlap;
815 
820  unsigned int n_sweeps;
821  };
822 
828  void initialize (const SparseMatrix &matrix,
829  const AdditionalData &additional_data = AdditionalData());
830  };
831 
832 
833 
870  {
871  public:
889  {
899  AdditionalData (const unsigned int ic_fill = 0,
900  const double ic_atol = 0.,
901  const double ic_rtol = 1.,
902  const unsigned int overlap = 0);
903 
912  unsigned int ic_fill;
913 
919  double ic_atol;
920 
925  double ic_rtol;
926 
931  unsigned int overlap;
932  };
933 
938  void initialize (const SparseMatrix &matrix,
939  const AdditionalData &additional_data = AdditionalData());
940  };
941 
942 
943 
974  {
975  public:
1014  {
1018  AdditionalData (const unsigned int ilu_fill = 0,
1019  const double ilu_atol = 0.,
1020  const double ilu_rtol = 1.,
1021  const unsigned int overlap = 0);
1022 
1026  unsigned int ilu_fill;
1027 
1032  double ilu_atol;
1033 
1038  double ilu_rtol;
1039 
1043  unsigned int overlap;
1044  };
1045 
1050  void initialize (const SparseMatrix &matrix,
1051  const AdditionalData &additional_data = AdditionalData());
1052  };
1053 
1054 
1055 
1056 
1057 
1058 
1095  {
1096  public:
1115  {
1126  AdditionalData (const double ilut_drop = 0.,
1127  const unsigned int ilut_fill = 0,
1128  const double ilut_atol = 0.,
1129  const double ilut_rtol = 1.,
1130  const unsigned int overlap = 0);
1131 
1136  double ilut_drop;
1137 
1146  unsigned int ilut_fill;
1147 
1153  double ilut_atol;
1154 
1159  double ilut_rtol;
1160 
1165  unsigned int overlap;
1166  };
1167 
1172  void initialize (const SparseMatrix &matrix,
1173  const AdditionalData &additional_data = AdditionalData());
1174  };
1175 
1176 
1177 
1197  {
1198  public:
1204  {
1208  AdditionalData (const unsigned int overlap = 0);
1209 
1210 
1215  unsigned int overlap;
1216  };
1217 
1222  void initialize (const SparseMatrix &matrix,
1223  const AdditionalData &additional_data = AdditionalData());
1224  };
1225 
1226 
1227 
1228 
1229 
1230 
1241  {
1242  public:
1248  {
1252  AdditionalData (const unsigned int degree = 1,
1253  const double max_eigenvalue = 10.,
1254  const double eigenvalue_ratio = 30.,
1255  const double min_eigenvalue = 1.,
1256  const double min_diagonal = 1e-12,
1257  const bool nonzero_starting = false);
1258 
1264  unsigned int degree;
1265 
1271 
1276 
1282 
1288 
1300  };
1301 
1306  void initialize (const SparseMatrix &matrix,
1307  const AdditionalData &additional_data = AdditionalData());
1308  };
1309 
1310 
1311 
1355  {
1356  public:
1357 
1365  {
1370  AdditionalData (const bool elliptic = true,
1371  const bool higher_order_elements = false,
1372  const unsigned int n_cycles = 1,
1373  const bool w_cyle = false,
1374  const double aggregation_threshold = 1e-4,
1375  const std::vector<std::vector<bool> > &constant_modes = std::vector<std::vector<bool> > (0),
1376  const unsigned int smoother_sweeps = 2,
1377  const unsigned int smoother_overlap = 0,
1378  const bool output_details = false,
1379  const char *smoother_type = "Chebyshev",
1380  const char *coarse_type = "Amesos-KLU");
1381 
1389  bool elliptic;
1390 
1396 
1401  unsigned int n_cycles;
1402 
1407  bool w_cycle;
1408 
1419 
1436  std::vector<std::vector<bool> > constant_modes;
1437 
1448  unsigned int smoother_sweeps;
1449 
1454  unsigned int smoother_overlap;
1455 
1462 
1497  const char *smoother_type;
1498 
1503  const char *coarse_type;
1504  };
1505 
1509  ~PreconditionAMG();
1510 
1511 
1517  void initialize (const SparseMatrix &matrix,
1518  const AdditionalData &additional_data = AdditionalData());
1519 
1538  void initialize (const Epetra_RowMatrix &matrix,
1539  const AdditionalData &additional_data = AdditionalData());
1540 
1553  void initialize (const SparseMatrix &matrix,
1554  const Teuchos::ParameterList &ml_parameters);
1555 
1563  void initialize (const Epetra_RowMatrix &matrix,
1564  const Teuchos::ParameterList &ml_parameters);
1565 
1572  template <typename number>
1573  void initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1574  const AdditionalData &additional_data = AdditionalData(),
1575  const double drop_tolerance = 1e-13,
1576  const ::SparsityPattern *use_this_sparsity = 0);
1577 
1590  void reinit ();
1591 
1596  void clear ();
1597 
1601  size_type memory_consumption () const;
1602 
1603  private:
1607  std_cxx11::shared_ptr<SparseMatrix> trilinos_matrix;
1608  };
1609 
1610 
1611 
1612 #if defined(DOXYGEN) || DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
1613 
1631  {
1632  public:
1633 
1634 
1642  {
1647  AdditionalData (const bool elliptic = true,
1648  const unsigned int n_cycles = 1,
1649  const bool w_cyle = false,
1650  const double aggregation_threshold = 1e-4,
1651  const std::vector<std::vector<bool> > &constant_modes = std::vector<std::vector<bool> > (0),
1652  const unsigned int smoother_sweeps = 2,
1653  const unsigned int smoother_overlap = 0,
1654  const bool output_details = false,
1655  const char *smoother_type = "Chebyshev",
1656  const char *coarse_type = "Amesos-KLU");
1657 
1665  bool elliptic;
1666 
1671  unsigned int n_cycles;
1672 
1677  bool w_cycle;
1678 
1689 
1696  std::vector<std::vector<bool> > constant_modes;
1697 
1708  unsigned int smoother_sweeps;
1709 
1714  unsigned int smoother_overlap;
1715 
1722 
1757  const char *smoother_type;
1758 
1763  const char *coarse_type;
1764  };
1765 
1770 
1775 
1781  void initialize (const SparseMatrix &matrix,
1782  const AdditionalData &additional_data = AdditionalData());
1783 
1790  void initialize (const Epetra_CrsMatrix &matrix,
1791  const AdditionalData &additional_data = AdditionalData());
1792 
1804  void initialize (const SparseMatrix &matrix,
1805  Teuchos::ParameterList &muelu_parameters);
1806 
1812  void initialize (const Epetra_CrsMatrix &matrix,
1813  Teuchos::ParameterList &muelu_parameters);
1814 
1821  template <typename number>
1822  void initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1823  const AdditionalData &additional_data = AdditionalData(),
1824  const double drop_tolerance = 1e-13,
1825  const ::SparsityPattern *use_this_sparsity = 0);
1826 
1831  void clear ();
1832 
1836  size_type memory_consumption () const;
1837 
1838  private:
1842  std_cxx11::shared_ptr<SparseMatrix> trilinos_matrix;
1843  };
1844 #endif
1845 
1846 
1847 
1857  {
1858  public:
1859 
1865  {
1870  };
1871 
1878  void initialize (const SparseMatrix &matrix,
1879  const AdditionalData &additional_data = AdditionalData());
1880 
1884  void vmult (VectorBase &dst,
1885  const VectorBase &src) const;
1886 
1890  void Tvmult (VectorBase &dst,
1891  const VectorBase &src) const;
1892 
1897  void vmult (::Vector<double> &dst,
1898  const ::Vector<double> &src) const;
1899 
1904  void Tvmult (::Vector<double> &dst,
1905  const ::Vector<double> &src) const;
1906 
1912  const ::LinearAlgebra::distributed::Vector<double> &src) const;
1913 
1920  const ::LinearAlgebra::distributed::Vector<double> &src) const;
1921  };
1922 
1923 
1924 
1925 // -------------------------- inline and template functions ----------------------
1926 
1927 
1928 #ifndef DOXYGEN
1929 
1930 
1931  inline
1932  void
1934  {
1935  // This only flips a flag that tells
1936  // Trilinos that any vmult operation
1937  // should be done with the
1938  // transpose. However, the matrix
1939  // structure is not reset.
1940  int ierr;
1941 
1942  if (!preconditioner->UseTranspose())
1943  {
1944  ierr = preconditioner->SetUseTranspose (true);
1945  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1946  }
1947  else
1948  {
1949  ierr = preconditioner->SetUseTranspose (false);
1950  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1951  }
1952  }
1953 
1954 
1955  inline
1956  void
1957  PreconditionBase::vmult (VectorBase &dst,
1958  const VectorBase &src) const
1959  {
1960  Assert (dst.vector_partitioner().SameAs(preconditioner->OperatorRangeMap()),
1961  ExcNonMatchingMaps("dst"));
1962  Assert (src.vector_partitioner().SameAs(preconditioner->OperatorDomainMap()),
1963  ExcNonMatchingMaps("src"));
1964 
1965  const int ierr = preconditioner->ApplyInverse (src.trilinos_vector(),
1966  dst.trilinos_vector());
1967  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1968  }
1969 
1970  inline
1971  void
1972  PreconditionBase::Tvmult (VectorBase &dst,
1973  const VectorBase &src) const
1974  {
1975  Assert (dst.vector_partitioner().SameAs(preconditioner->OperatorRangeMap()),
1976  ExcNonMatchingMaps("dst"));
1977  Assert (src.vector_partitioner().SameAs(preconditioner->OperatorDomainMap()),
1978  ExcNonMatchingMaps("src"));
1979 
1980  preconditioner->SetUseTranspose(true);
1981  const int ierr = preconditioner->ApplyInverse (src.trilinos_vector(),
1982  dst.trilinos_vector());
1983  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1984  preconditioner->SetUseTranspose(false);
1985  }
1986 
1987 
1988  // For the implementation of the <code>vmult</code> function with deal.II
1989  // data structures we note that invoking a call of the Trilinos
1990  // preconditioner requires us to use Epetra vectors as well. We do this by
1991  // providing a view, i.e., feed Trilinos with a pointer to the data, so we
1992  // avoid copying the content of the vectors during the iteration (this
1993  // function is only useful when used in serial anyway). In the declaration
1994  // of the right hand side, we need to cast the source vector (that is
1995  // <code>const</code> in all deal.II calls) to non-constant value, as this
1996  // is the way Trilinos wants to have them.
1997  inline
1999  const ::Vector<double> &src) const
2000  {
2001  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.size()),
2002  preconditioner->OperatorDomainMap().NumMyElements());
2003  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.size()),
2004  preconditioner->OperatorRangeMap().NumMyElements());
2005  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
2006  dst.begin());
2007  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
2008  const_cast<double *>(src.begin()));
2009 
2010  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
2011  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
2012  }
2013 
2014 
2015  inline
2017  const ::Vector<double> &src) const
2018  {
2019  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.size()),
2020  preconditioner->OperatorDomainMap().NumMyElements());
2021  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.size()),
2022  preconditioner->OperatorRangeMap().NumMyElements());
2023  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
2024  dst.begin());
2025  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
2026  const_cast<double *>(src.begin()));
2027 
2028  preconditioner->SetUseTranspose(true);
2029  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
2030  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
2031  preconditioner->SetUseTranspose(false);
2032  }
2033 
2034 
2035 
2036  inline
2037  void
2040  {
2041  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.local_size()),
2042  preconditioner->OperatorDomainMap().NumMyElements());
2043  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.local_size()),
2044  preconditioner->OperatorRangeMap().NumMyElements());
2045  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
2046  dst.begin());
2047  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
2048  const_cast<double *>(src.begin()));
2049 
2050  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
2051  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
2052  }
2053 
2054  inline
2055  void
2058  {
2059  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.local_size()),
2060  preconditioner->OperatorDomainMap().NumMyElements());
2061  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.local_size()),
2062  preconditioner->OperatorRangeMap().NumMyElements());
2063  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
2064  dst.begin());
2065  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
2066  const_cast<double *>(src.begin()));
2067 
2068  preconditioner->SetUseTranspose(true);
2069  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
2070  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
2071  preconditioner->SetUseTranspose(false);
2072  }
2073 
2074 #endif
2075 
2076 }
2077 
2078 
2082 DEAL_II_NAMESPACE_CLOSE
2083 
2084 #endif // DEAL_II_WITH_TRILINOS
2085 
2086 /*---------------------------- trilinos_precondition.h ---------------------------*/
2087 
2088 #endif
2089 /*---------------------------- trilinos_precondition.h ---------------------------*/
static ::ExceptionBase & ExcTrilinosError(int arg1)
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const bool elliptic=true, const unsigned int n_cycles=1, const bool w_cyle=false, const double aggregation_threshold=1e-4, const std::vector< std::vector< bool > > &constant_modes=std::vector< std::vector< bool > >(0), const unsigned int smoother_sweeps=2, const unsigned int smoother_overlap=0, const bool output_details=false, const char *smoother_type="Chebyshev", const char *coarse_type="Amesos-KLU")
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
virtual void Tvmult(VectorBase &dst, const VectorBase &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
AdditionalData(const unsigned int ilu_fill=0, const double ilu_atol=0., const double ilu_rtol=1., const unsigned int overlap=0)
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
virtual void vmult(VectorBase &dst, const VectorBase &src) const
std_cxx11::shared_ptr< Epetra_Map > vector_distributor
std_cxx11::shared_ptr< Epetra_Operator > preconditioner
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
#define Assert(cond, exc)
Definition: exceptions.h:313
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
AdditionalData(const unsigned int degree=1, const double max_eigenvalue=10., const double eigenvalue_ratio=30., const double min_eigenvalue=1., const double min_diagonal=1e-12, const bool nonzero_starting=false)
AdditionalData(const bool elliptic=true, const bool higher_order_elements=false, const unsigned int n_cycles=1, const bool w_cyle=false, const double aggregation_threshold=1e-4, const std::vector< std::vector< bool > > &constant_modes=std::vector< std::vector< bool > >(0), const unsigned int smoother_sweeps=2, const unsigned int smoother_overlap=0, const bool output_details=false, const char *smoother_type="Chebyshev", const char *coarse_type="Amesos-KLU")
std::size_t size() const
AdditionalData(const double ilut_drop=0., const unsigned int ilut_fill=0, const double ilut_atol=0., const double ilut_rtol=1., const unsigned int overlap=0)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
iterator begin()
Epetra_Operator & trilinos_operator() const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void vmult(VectorBase &dst, const VectorBase &src) const
void Tvmult(VectorBase &dst, const VectorBase &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
static ::ExceptionBase & ExcNonMatchingMaps(std::string arg1)
AdditionalData(const unsigned int ic_fill=0, const double ic_atol=0., const double ic_rtol=1., const unsigned int overlap=0)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const unsigned int block_size=1, const std::string &block_creation_type="linear", const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)