Reference documentation for deal.II version 8.4.2
trilinos_precondition.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2008 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__trilinos_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/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>
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 
78  {
79  public:
83  typedef ::types::global_dof_index size_type;
84 
90  {};
91 
98 
103 
107  ~PreconditionBase ();
108 
113  void clear ();
114 
118  virtual void vmult (VectorBase &dst,
119  const VectorBase &src) const;
120 
124  virtual void Tvmult (VectorBase &dst,
125  const VectorBase &src) const;
126 
131  virtual void vmult (::Vector<double> &dst,
132  const ::Vector<double> &src) const;
133 
138  virtual void Tvmult (::Vector<double> &dst,
139  const ::Vector<double> &src) const;
140 
145  virtual void vmult (::parallel::distributed::Vector<double> &dst,
146  const ::parallel::distributed::Vector<double> &src) const;
147 
152  virtual void Tvmult (::parallel::distributed::Vector<double> &dst,
153  const ::parallel::distributed::Vector<double> &src) const;
154 
162  Epetra_Operator &trilinos_operator() const;
163 
167  DeclException1 (ExcNonMatchingMaps,
168  std::string,
169  << "The sparse matrix the preconditioner is based on "
170  << "uses a map that is not compatible to the one in vector "
171  << arg1
172  << ". Check preconditioner and matrix setup.");
173 
174  friend class SolverBase;
175 
176  protected:
181  std_cxx11::shared_ptr<Epetra_Operator> preconditioner;
182 
187 #ifdef DEAL_II_WITH_MPI
188  Epetra_MpiComm communicator;
189 #else
190  Epetra_SerialComm communicator;
191 #endif
192 
197  std_cxx11::shared_ptr<Epetra_Map> vector_distributor;
198  };
199 
200 
218  {
219  public:
220 
233  {
238  AdditionalData (const double omega = 1,
239  const double min_diagonal = 0,
240  const unsigned int n_sweeps = 1);
241 
245  double omega;
246 
254  double min_diagonal;
255 
260  unsigned int n_sweeps;
261  };
262 
267  void initialize (const SparseMatrix &matrix,
268  const AdditionalData &additional_data = AdditionalData());
269  };
270 
271 
272 
273 
301  {
302  public:
303 
318  {
325  AdditionalData (const double omega = 1,
326  const double min_diagonal = 0,
327  const unsigned int overlap = 0,
328  const unsigned int n_sweeps = 1);
329 
334  double omega;
335 
343  double min_diagonal;
344 
349  unsigned int overlap;
350 
355  unsigned int n_sweeps;
356  };
357 
363  void initialize (const SparseMatrix &matrix,
364  const AdditionalData &additional_data = AdditionalData());
365  };
366 
367 
368 
369 
397  {
398  public:
399 
414  {
421  AdditionalData (const double omega = 1,
422  const double min_diagonal = 0,
423  const unsigned int overlap = 0,
424  const unsigned int n_sweeps = 1);
425 
430  double omega;
431 
439  double min_diagonal;
440 
445  unsigned int overlap;
446 
451  unsigned int n_sweeps;
452  };
453 
459  void initialize (const SparseMatrix &matrix,
460  const AdditionalData &additional_data = AdditionalData());
461  };
462 
463 
464 
483  {
484  public:
485 
504  {
510  AdditionalData (const unsigned int block_size = 1,
511  const std::string block_creation_type = "linear",
512  const double omega = 1,
513  const double min_diagonal = 0,
514  const unsigned int n_sweeps = 1);
515 
519  unsigned int block_size;
520 
528  std::string block_creation_type;
529 
534  double omega;
535 
543  double min_diagonal;
544 
549  unsigned int n_sweeps;
550  };
551 
556  void initialize (const SparseMatrix &matrix,
557  const AdditionalData &additional_data = AdditionalData());
558  };
559 
560 
561 
562 
586  {
587  public:
588 
607  {
615  AdditionalData (const unsigned int block_size = 1,
616  const std::string block_creation_type = "linear",
617  const double omega = 1,
618  const double min_diagonal = 0,
619  const unsigned int overlap = 0,
620  const unsigned int n_sweeps = 1);
621 
625  unsigned int block_size;
626 
634  std::string block_creation_type;
635 
640  double omega;
641 
649  double min_diagonal;
650 
655  unsigned int overlap;
656 
661  unsigned int n_sweeps;
662  };
663 
669  void initialize (const SparseMatrix &matrix,
670  const AdditionalData &additional_data = AdditionalData());
671  };
672 
673 
674 
675 
699  {
700  public:
701 
720  {
728  AdditionalData (const unsigned int block_size = 1,
729  const std::string block_creation_type = "linear",
730  const double omega = 1,
731  const double min_diagonal = 0,
732  const unsigned int overlap = 0,
733  const unsigned int n_sweeps = 1);
734 
738  unsigned int block_size;
739 
747  std::string block_creation_type;
748 
753  double omega;
754 
762  double min_diagonal;
763 
768  unsigned int overlap;
769 
774  unsigned int n_sweeps;
775  };
776 
782  void initialize (const SparseMatrix &matrix,
783  const AdditionalData &additional_data = AdditionalData());
784  };
785 
786 
787 
824  {
825  public:
843  {
853  AdditionalData (const unsigned int ic_fill = 0,
854  const double ic_atol = 0.,
855  const double ic_rtol = 1.,
856  const unsigned int overlap = 0);
857 
866  unsigned int ic_fill;
867 
873  double ic_atol;
874 
879  double ic_rtol;
880 
885  unsigned int overlap;
886  };
887 
892  void initialize (const SparseMatrix &matrix,
893  const AdditionalData &additional_data = AdditionalData());
894  };
895 
896 
897 
928  {
929  public:
968  {
972  AdditionalData (const unsigned int ilu_fill = 0,
973  const double ilu_atol = 0.,
974  const double ilu_rtol = 1.,
975  const unsigned int overlap = 0);
976 
980  unsigned int ilu_fill;
981 
986  double ilu_atol;
987 
992  double ilu_rtol;
993 
997  unsigned int overlap;
998  };
999 
1004  void initialize (const SparseMatrix &matrix,
1005  const AdditionalData &additional_data = AdditionalData());
1006  };
1007 
1008 
1009 
1010 
1011 
1012 
1049  {
1050  public:
1069  {
1080  AdditionalData (const double ilut_drop = 0.,
1081  const unsigned int ilut_fill = 0,
1082  const double ilut_atol = 0.,
1083  const double ilut_rtol = 1.,
1084  const unsigned int overlap = 0);
1085 
1090  double ilut_drop;
1091 
1100  unsigned int ilut_fill;
1101 
1107  double ilut_atol;
1108 
1113  double ilut_rtol;
1114 
1119  unsigned int overlap;
1120  };
1121 
1126  void initialize (const SparseMatrix &matrix,
1127  const AdditionalData &additional_data = AdditionalData());
1128  };
1129 
1130 
1131 
1151  {
1152  public:
1158  {
1162  AdditionalData (const unsigned int overlap = 0);
1163 
1164 
1169  unsigned int overlap;
1170  };
1171 
1176  void initialize (const SparseMatrix &matrix,
1177  const AdditionalData &additional_data = AdditionalData());
1178  };
1179 
1180 
1181 
1182 
1183 
1184 
1195  {
1196  public:
1202  {
1206  AdditionalData (const unsigned int degree = 1,
1207  const double max_eigenvalue = 10.,
1208  const double eigenvalue_ratio = 30.,
1209  const double min_eigenvalue = 1.,
1210  const double min_diagonal = 1e-12,
1211  const bool nonzero_starting = false);
1212 
1218  unsigned int degree;
1219 
1225 
1230 
1236 
1242 
1254  };
1255 
1260  void initialize (const SparseMatrix &matrix,
1261  const AdditionalData &additional_data = AdditionalData());
1262  };
1263 
1264 
1265 
1309  {
1310  public:
1311 
1319  {
1324  AdditionalData (const bool elliptic = true,
1325  const bool higher_order_elements = false,
1326  const unsigned int n_cycles = 1,
1327  const bool w_cyle = false,
1328  const double aggregation_threshold = 1e-4,
1329  const std::vector<std::vector<bool> > &constant_modes = std::vector<std::vector<bool> > (0),
1330  const unsigned int smoother_sweeps = 2,
1331  const unsigned int smoother_overlap = 0,
1332  const bool output_details = false,
1333  const char *smoother_type = "Chebyshev",
1334  const char *coarse_type = "Amesos-KLU");
1335 
1343  bool elliptic;
1344 
1350 
1355  unsigned int n_cycles;
1356 
1361  bool w_cycle;
1362 
1373 
1390  std::vector<std::vector<bool> > constant_modes;
1391 
1402  unsigned int smoother_sweeps;
1403 
1408  unsigned int smoother_overlap;
1409 
1416 
1451  const char *smoother_type;
1452 
1457  const char *coarse_type;
1458  };
1459 
1463  ~PreconditionAMG();
1464 
1465 
1471  void initialize (const SparseMatrix &matrix,
1472  const AdditionalData &additional_data = AdditionalData());
1473 
1492  void initialize (const Epetra_RowMatrix &matrix,
1493  const AdditionalData &additional_data = AdditionalData());
1494 
1507  void initialize (const SparseMatrix &matrix,
1508  const Teuchos::ParameterList &ml_parameters);
1509 
1517  void initialize (const Epetra_RowMatrix &matrix,
1518  const Teuchos::ParameterList &ml_parameters);
1519 
1526  template <typename number>
1527  void initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1528  const AdditionalData &additional_data = AdditionalData(),
1529  const double drop_tolerance = 1e-13,
1530  const ::SparsityPattern *use_this_sparsity = 0);
1531 
1544  void reinit ();
1545 
1550  void clear ();
1551 
1555  size_type memory_consumption () const;
1556 
1557  private:
1561  std_cxx11::shared_ptr<SparseMatrix> trilinos_matrix;
1562  };
1563 
1564 
1565 
1566 #if defined(DOXYGEN) || DEAL_II_TRILINOS_VERSION_GTE(11,14,0)
1567 
1584  {
1585  public:
1586 
1587 
1595  {
1600  AdditionalData (const bool elliptic = true,
1601  const unsigned int n_cycles = 1,
1602  const bool w_cyle = false,
1603  const double aggregation_threshold = 1e-4,
1604  const std::vector<std::vector<bool> > &constant_modes = std::vector<std::vector<bool> > (0),
1605  const unsigned int smoother_sweeps = 2,
1606  const unsigned int smoother_overlap = 0,
1607  const bool output_details = false,
1608  const char *smoother_type = "Chebyshev",
1609  const char *coarse_type = "Amesos-KLU");
1610 
1618  bool elliptic;
1619 
1624  unsigned int n_cycles;
1625 
1630  bool w_cycle;
1631 
1642 
1649  std::vector<std::vector<bool> > constant_modes;
1650 
1661  unsigned int smoother_sweeps;
1662 
1667  unsigned int smoother_overlap;
1668 
1675 
1710  const char *smoother_type;
1711 
1716  const char *coarse_type;
1717  };
1718 
1723 
1729  void initialize (const SparseMatrix &matrix,
1730  const AdditionalData &additional_data = AdditionalData());
1731 
1738  void initialize (const Epetra_CrsMatrix &matrix,
1739  const AdditionalData &additional_data = AdditionalData());
1740 
1752  void initialize (const SparseMatrix &matrix,
1753  Teuchos::ParameterList &muelu_parameters);
1754 
1760  void initialize (const Epetra_CrsMatrix &matrix,
1761  Teuchos::ParameterList &muelu_parameters);
1762 
1769  template <typename number>
1770  void initialize (const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1771  const AdditionalData &additional_data = AdditionalData(),
1772  const double drop_tolerance = 1e-13,
1773  const ::SparsityPattern *use_this_sparsity = 0);
1774 
1779  void clear ();
1780 
1784  size_type memory_consumption () const;
1785 
1786  private:
1790  std_cxx11::shared_ptr<SparseMatrix> trilinos_matrix;
1791  };
1792 #endif
1793 
1794 
1795 
1804  {
1805  public:
1806 
1810  void vmult (VectorBase &dst,
1811  const VectorBase &src) const;
1812 
1816  void Tvmult (VectorBase &dst,
1817  const VectorBase &src) const;
1818 
1823  void vmult (::Vector<double> &dst,
1824  const ::Vector<double> &src) const;
1825 
1830  void Tvmult (::Vector<double> &dst,
1831  const ::Vector<double> &src) const;
1832 
1837  void vmult (parallel::distributed::Vector<double> &dst,
1838  const ::parallel::distributed::Vector<double> &src) const;
1839 
1845  void Tvmult (parallel::distributed::Vector<double> &dst,
1846  const ::parallel::distributed::Vector<double> &src) const;
1847  };
1848 
1849 
1850 
1851 // -------------------------- inline and template functions ----------------------
1852 
1853 
1854 #ifndef DOXYGEN
1855 
1856  inline
1857  void
1858  PreconditionBase::vmult (VectorBase &dst,
1859  const VectorBase &src) const
1860  {
1861  Assert (dst.vector_partitioner().SameAs(preconditioner->OperatorRangeMap()),
1862  ExcNonMatchingMaps("dst"));
1863  Assert (src.vector_partitioner().SameAs(preconditioner->OperatorDomainMap()),
1864  ExcNonMatchingMaps("src"));
1865 
1866  const int ierr = preconditioner->ApplyInverse (src.trilinos_vector(),
1867  dst.trilinos_vector());
1868  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1869  }
1870 
1871  inline
1872  void
1873  PreconditionBase::Tvmult (VectorBase &dst,
1874  const VectorBase &src) const
1875  {
1876  Assert (dst.vector_partitioner().SameAs(preconditioner->OperatorRangeMap()),
1877  ExcNonMatchingMaps("dst"));
1878  Assert (src.vector_partitioner().SameAs(preconditioner->OperatorDomainMap()),
1879  ExcNonMatchingMaps("src"));
1880 
1881  preconditioner->SetUseTranspose(true);
1882  const int ierr = preconditioner->ApplyInverse (src.trilinos_vector(),
1883  dst.trilinos_vector());
1884  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1885  preconditioner->SetUseTranspose(false);
1886  }
1887 
1888 
1889  // For the implementation of the <code>vmult</code> function with deal.II
1890  // data structures we note that invoking a call of the Trilinos
1891  // preconditioner requires us to use Epetra vectors as well. We do this by
1892  // providing a view, i.e., feed Trilinos with a pointer to the data, so we
1893  // avoid copying the content of the vectors during the iteration (this
1894  // function is only useful when used in serial anyway). In the declaration
1895  // of the right hand side, we need to cast the source vector (that is
1896  // <code>const</code> in all deal.II calls) to non-constant value, as this
1897  // is the way Trilinos wants to have them.
1898  inline
1899  void PreconditionBase::vmult (::Vector<double> &dst,
1900  const ::Vector<double> &src) const
1901  {
1902  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.size()),
1903  preconditioner->OperatorDomainMap().NumMyElements());
1904  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.size()),
1905  preconditioner->OperatorRangeMap().NumMyElements());
1906  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
1907  dst.begin());
1908  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
1909  const_cast<double *>(src.begin()));
1910 
1911  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
1912  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1913  }
1914 
1915 
1916  inline
1917  void PreconditionBase::Tvmult (::Vector<double> &dst,
1918  const ::Vector<double> &src) const
1919  {
1920  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.size()),
1921  preconditioner->OperatorDomainMap().NumMyElements());
1922  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.size()),
1923  preconditioner->OperatorRangeMap().NumMyElements());
1924  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
1925  dst.begin());
1926  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
1927  const_cast<double *>(src.begin()));
1928 
1929  preconditioner->SetUseTranspose(true);
1930  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
1931  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1932  preconditioner->SetUseTranspose(false);
1933  }
1934 
1935 
1936 
1937  inline
1938  void
1939  PreconditionBase::vmult (parallel::distributed::Vector<double> &dst,
1940  const parallel::distributed::Vector<double> &src) const
1941  {
1942  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.local_size()),
1943  preconditioner->OperatorDomainMap().NumMyElements());
1944  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.local_size()),
1945  preconditioner->OperatorRangeMap().NumMyElements());
1946  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
1947  dst.begin());
1948  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
1949  const_cast<double *>(src.begin()));
1950 
1951  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
1952  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1953  }
1954 
1955  inline
1956  void
1957  PreconditionBase::Tvmult (parallel::distributed::Vector<double> &dst,
1958  const parallel::distributed::Vector<double> &src) const
1959  {
1960  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(dst.local_size()),
1961  preconditioner->OperatorDomainMap().NumMyElements());
1962  AssertDimension (static_cast<TrilinosWrappers::types::int_type>(src.local_size()),
1963  preconditioner->OperatorRangeMap().NumMyElements());
1964  Epetra_Vector tril_dst (View, preconditioner->OperatorDomainMap(),
1965  dst.begin());
1966  Epetra_Vector tril_src (View, preconditioner->OperatorRangeMap(),
1967  const_cast<double *>(src.begin()));
1968 
1969  preconditioner->SetUseTranspose(true);
1970  const int ierr = preconditioner->ApplyInverse (tril_src, tril_dst);
1971  AssertThrow (ierr == 0, ExcTrilinosError(ierr));
1972  preconditioner->SetUseTranspose(false);
1973  }
1974 
1975  inline
1976  Epetra_Operator &
1977  PreconditionBase::trilinos_operator () const
1978  {
1979  AssertThrow (preconditioner, ExcMessage("Trying to dereference a null pointer."));
1980  return (*preconditioner);
1981  }
1982 
1983 #endif
1984 
1985 }
1986 
1987 
1991 DEAL_II_NAMESPACE_CLOSE
1992 
1993 #endif // DEAL_II_WITH_TRILINOS
1994 
1995 /*---------------------------- trilinos_precondition.h ---------------------------*/
1996 
1997 #endif
1998 /*---------------------------- trilinos_precondition.h ---------------------------*/
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
std_cxx11::shared_ptr< Epetra_Map > vector_distributor
std_cxx11::shared_ptr< Epetra_Operator > preconditioner
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:542
#define Assert(cond, exc)
Definition: exceptions.h:294
std_cxx11::shared_ptr< SparseMatrix > trilinos_matrix
const Epetra_Map & vector_partitioner() const
std::size_t size() const
iterator begin()
const Epetra_MultiVector & trilinos_vector() const
size_type local_size() const