Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
trilinos_precondition.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2008 - 2022 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.md at
12// the top level directory of deal.II.
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
25
28
29# include <Epetra_Map.h>
30# include <Epetra_MpiComm.h>
31# include <Epetra_MultiVector.h>
32# include <Epetra_RowMatrix.h>
33# include <Epetra_Vector.h>
34# include <Teuchos_ParameterList.hpp>
35
36# include <memory>
37
38// forward declarations
39# ifndef DOXYGEN
40class Ifpack_Preconditioner;
41class Ifpack_Chebyshev;
42namespace ML_Epetra
43{
44 class MultiLevelPreconditioner;
45}
46# endif
47
49
50// forward declarations
51# ifndef DOXYGEN
52template <typename number>
53class SparseMatrix;
54template <typename number>
55class Vector;
56class SparsityPattern;
57# endif
58
63namespace TrilinosWrappers
64{
65 // forward declarations
66 class SparseMatrix;
68 class SolverBase;
69
77 {
78 public:
83
89 {};
90
97
102
106 ~PreconditionBase() override = default;
107
112 void
113 clear();
114
119 get_mpi_communicator() const;
120
130 void
132
136 virtual void
137 vmult(MPI::Vector &dst, const MPI::Vector &src) const;
138
142 virtual void
143 Tvmult(MPI::Vector &dst, const MPI::Vector &src) const;
144
149 virtual void
150 vmult(::Vector<double> &dst, const ::Vector<double> &src) const;
151
156 virtual void
158 const ::Vector<double> &src) const;
159
164 virtual void
166 const ::LinearAlgebra::distributed::Vector<double> &src) const;
167
172 virtual void
174 const ::LinearAlgebra::distributed::Vector<double> &src) const;
175
186 trilinos_operator() const;
188
193
200
208
210
219 std::string,
220 << "The sparse matrix the preconditioner is based on "
221 << "uses a map that is not compatible to the one in vector "
222 << arg1 << ". Check preconditioner and matrix setup.");
224
225 friend class SolverBase;
226
227 protected:
232 Teuchos::RCP<Epetra_Operator> preconditioner;
233
238 Epetra_MpiComm communicator;
239
244 std::shared_ptr<Epetra_Map> vector_distributor;
245 };
246
247
264 {
265 public:
278 {
283 AdditionalData(const double omega = 1,
284 const double min_diagonal = 0,
285 const unsigned int n_sweeps = 1);
286
290 double omega;
291
300
305 unsigned int n_sweeps;
306 };
307
312 void
313 initialize(const SparseMatrix & matrix,
314 const AdditionalData &additional_data = AdditionalData());
315 };
316
317
318
345 {
346 public:
361 {
368 AdditionalData(const double omega = 1,
369 const double min_diagonal = 0,
370 const unsigned int overlap = 0,
371 const unsigned int n_sweeps = 1);
372
377 double omega;
378
387
392 unsigned int overlap;
393
398 unsigned int n_sweeps;
399 };
400
406 void
407 initialize(const SparseMatrix & matrix,
408 const AdditionalData &additional_data = AdditionalData());
409 };
410
411
412
439 {
440 public:
455 {
462 AdditionalData(const double omega = 1,
463 const double min_diagonal = 0,
464 const unsigned int overlap = 0,
465 const unsigned int n_sweeps = 1);
466
471 double omega;
472
481
486 unsigned int overlap;
487
492 unsigned int n_sweeps;
493 };
494
500 void
501 initialize(const SparseMatrix & matrix,
502 const AdditionalData &additional_data = AdditionalData());
503 };
504
505
506
524 {
525 public:
542 {
548 AdditionalData(const unsigned int block_size = 1,
549 const std::string &block_creation_type = "linear",
550 const double omega = 1,
551 const double min_diagonal = 0,
552 const unsigned int n_sweeps = 1);
553
557 unsigned int block_size;
558
567
572 double omega;
573
582
587 unsigned int n_sweeps;
588 };
589
594 void
595 initialize(const SparseMatrix & matrix,
596 const AdditionalData &additional_data = AdditionalData());
597 };
598
599
600
623 {
624 public:
643 {
651 AdditionalData(const unsigned int block_size = 1,
652 const std::string &block_creation_type = "linear",
653 const double omega = 1,
654 const double min_diagonal = 0,
655 const unsigned int overlap = 0,
656 const unsigned int n_sweeps = 1);
657
661 unsigned int block_size;
662
671
676 double omega;
677
686
691 unsigned int overlap;
692
697 unsigned int n_sweeps;
698 };
699
705 void
706 initialize(const SparseMatrix & matrix,
707 const AdditionalData &additional_data = AdditionalData());
708 };
709
710
711
734 {
735 public:
754 {
762 AdditionalData(const unsigned int block_size = 1,
763 const std::string &block_creation_type = "linear",
764 const double omega = 1,
765 const double min_diagonal = 0,
766 const unsigned int overlap = 0,
767 const unsigned int n_sweeps = 1);
768
772 unsigned int block_size;
773
782
787 double omega;
788
797
802 unsigned int overlap;
803
808 unsigned int n_sweeps;
809 };
810
816 void
817 initialize(const SparseMatrix & matrix,
818 const AdditionalData &additional_data = AdditionalData());
819 };
820
821
822
858 {
859 public:
877 {
887 AdditionalData(const unsigned int ic_fill = 0,
888 const double ic_atol = 0.,
889 const double ic_rtol = 1.,
890 const unsigned int overlap = 0);
891
900 unsigned int ic_fill;
901
907 double ic_atol;
908
913 double ic_rtol;
914
919 unsigned int overlap;
920 };
921
926 void
927 initialize(const SparseMatrix & matrix,
928 const AdditionalData &additional_data = AdditionalData());
929 };
930
931
932
962 {
963 public:
1002 {
1006 AdditionalData(const unsigned int ilu_fill = 0,
1007 const double ilu_atol = 0.,
1008 const double ilu_rtol = 1.,
1009 const unsigned int overlap = 0);
1010
1014 unsigned int ilu_fill;
1015
1020 double ilu_atol;
1021
1026 double ilu_rtol;
1027
1031 unsigned int overlap;
1032 };
1033
1038 void
1039 initialize(const SparseMatrix & matrix,
1040 const AdditionalData &additional_data = AdditionalData());
1041 };
1042
1043
1044
1080 {
1081 public:
1100 {
1111 AdditionalData(const double ilut_drop = 0.,
1112 const unsigned int ilut_fill = 0,
1113 const double ilut_atol = 0.,
1114 const double ilut_rtol = 1.,
1115 const unsigned int overlap = 0);
1116
1122
1131 unsigned int ilut_fill;
1132
1139
1145
1150 unsigned int overlap;
1151 };
1152
1157 void
1158 initialize(const SparseMatrix & matrix,
1159 const AdditionalData &additional_data = AdditionalData());
1160 };
1161
1162
1163
1182 {
1183 public:
1189 {
1193 AdditionalData(const unsigned int overlap = 0);
1194
1195
1200 unsigned int overlap;
1201 };
1202
1207 void
1208 initialize(const SparseMatrix & matrix,
1209 const AdditionalData &additional_data = AdditionalData());
1210 };
1211
1212
1213
1223 {
1224 public:
1230 {
1234 AdditionalData(const unsigned int degree = 1,
1235 const double max_eigenvalue = 10.,
1236 const double eigenvalue_ratio = 30.,
1237 const double min_eigenvalue = 1.,
1238 const double min_diagonal = 1e-12,
1239 const bool nonzero_starting = false);
1240
1246 unsigned int degree;
1247
1253
1258
1264
1270
1282 };
1283
1288 void
1289 initialize(const SparseMatrix & matrix,
1290 const AdditionalData &additional_data = AdditionalData());
1291 };
1292
1293
1294
1337 {
1338 public:
1346 {
1370 AdditionalData(const bool elliptic = true,
1371 const bool higher_order_elements = false,
1372 const unsigned int n_cycles = 1,
1373 const bool w_cycle = false,
1374 const double aggregation_threshold = 1e-4,
1375 const std::vector<std::vector<bool>> &constant_modes =
1376 std::vector<std::vector<bool>>(0),
1377 const unsigned int smoother_sweeps = 2,
1378 const unsigned int smoother_overlap = 0,
1379 const bool output_details = false,
1380 const char * smoother_type = "Chebyshev",
1381 const char * coarse_type = "Amesos-KLU");
1382
1413 void
1415 Teuchos::ParameterList & parameter_list,
1416 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1417 const Epetra_RowMatrix & matrix) const;
1418
1426 void
1428 Teuchos::ParameterList & parameter_list,
1429 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1430 const SparseMatrix & matrix) const;
1431
1437 void
1439 Teuchos::ParameterList & parameter_list,
1440 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1441 const Epetra_RowMatrix & matrix) const;
1442
1448 void
1450 Teuchos::ParameterList & parameter_list,
1451 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1452 const SparseMatrix & matrix) const;
1453
1462
1468
1473 unsigned int n_cycles;
1474
1480
1491
1516 std::vector<std::vector<bool>> constant_modes;
1517
1528 unsigned int smoother_sweeps;
1529
1534 unsigned int smoother_overlap;
1535
1542
1577 const char *smoother_type;
1578
1583 const char *coarse_type;
1584 };
1585
1589 ~PreconditionAMG() override;
1590
1591
1597 void
1598 initialize(const SparseMatrix & matrix,
1599 const AdditionalData &additional_data = AdditionalData());
1600
1619 void
1620 initialize(const Epetra_RowMatrix &matrix,
1621 const AdditionalData & additional_data = AdditionalData());
1622
1637 void
1638 initialize(const SparseMatrix & matrix,
1639 const Teuchos::ParameterList &ml_parameters);
1640
1648 void
1649 initialize(const Epetra_RowMatrix & matrix,
1650 const Teuchos::ParameterList &ml_parameters);
1651
1658 template <typename number>
1659 void
1660 initialize(const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1661 const AdditionalData &additional_data = AdditionalData(),
1662 const double drop_tolerance = 1e-13,
1663 const ::SparsityPattern *use_this_sparsity = nullptr);
1664
1677 void
1678 reinit();
1679
1684 void
1685 clear();
1686
1690 size_type
1691 memory_consumption() const;
1692
1693 private:
1697 std::shared_ptr<SparseMatrix> trilinos_matrix;
1698 };
1699
1700
1701
1702# if defined(DOXYGEN) || defined(DEAL_II_TRILINOS_WITH_MUELU)
1722 {
1723 public:
1731 {
1736 AdditionalData(const bool elliptic = true,
1737 const unsigned int n_cycles = 1,
1738 const bool w_cycle = false,
1739 const double aggregation_threshold = 1e-4,
1740 const std::vector<std::vector<bool>> &constant_modes =
1741 std::vector<std::vector<bool>>(0),
1742 const unsigned int smoother_sweeps = 2,
1743 const unsigned int smoother_overlap = 0,
1744 const bool output_details = false,
1745 const char * smoother_type = "Chebyshev",
1746 const char * coarse_type = "Amesos-KLU");
1747
1756
1761 unsigned int n_cycles;
1762
1768
1779
1786 std::vector<std::vector<bool>> constant_modes;
1787
1798 unsigned int smoother_sweeps;
1799
1804 unsigned int smoother_overlap;
1805
1812
1847 const char *smoother_type;
1848
1853 const char *coarse_type;
1854 };
1855
1860
1864 virtual ~PreconditionAMGMueLu() override = default;
1865
1871 void
1872 initialize(const SparseMatrix & matrix,
1873 const AdditionalData &additional_data = AdditionalData());
1874
1881 void
1882 initialize(const Epetra_CrsMatrix &matrix,
1883 const AdditionalData & additional_data = AdditionalData());
1884
1898 void
1899 initialize(const SparseMatrix & matrix,
1900 Teuchos::ParameterList &muelu_parameters);
1901
1907 void
1908 initialize(const Epetra_CrsMatrix &matrix,
1909 Teuchos::ParameterList &muelu_parameters);
1910
1917 template <typename number>
1918 void
1919 initialize(const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1920 const AdditionalData &additional_data = AdditionalData(),
1921 const double drop_tolerance = 1e-13,
1922 const ::SparsityPattern *use_this_sparsity = nullptr);
1923
1928 void
1929 clear();
1930
1934 size_type
1935 memory_consumption() const;
1936
1937 private:
1941 std::shared_ptr<SparseMatrix> trilinos_matrix;
1942 };
1943# endif
1944
1945
1946
1954 {
1955 public:
1961 {};
1962
1969 void
1970 initialize(const SparseMatrix & matrix,
1971 const AdditionalData &additional_data = AdditionalData());
1972
1976 void
1977 vmult(MPI::Vector &dst, const MPI::Vector &src) const override;
1978
1982 void
1983 Tvmult(MPI::Vector &dst, const MPI::Vector &src) const override;
1984
1989 void
1990 vmult(::Vector<double> & dst,
1991 const ::Vector<double> &src) const override;
1992
1997 void
1999 const ::Vector<double> &src) const override;
2000
2005 void
2007 const ::LinearAlgebra::distributed::Vector<double> &src)
2008 const override;
2009
2015 void
2017 const ::LinearAlgebra::distributed::Vector<double> &src)
2018 const override;
2019 };
2020
2021
2022
2023 // ----------------------- inline and template functions --------------------
2024
2025
2026# ifndef DOXYGEN
2027
2028
2029 inline void
2031 {
2032 // This only flips a flag that tells
2033 // Trilinos that any vmult operation
2034 // should be done with the
2035 // transpose. However, the matrix
2036 // structure is not reset.
2037 int ierr;
2038
2039 if (!preconditioner->UseTranspose())
2040 {
2041 ierr = preconditioner->SetUseTranspose(true);
2042 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2043 }
2044 else
2045 {
2046 ierr = preconditioner->SetUseTranspose(false);
2047 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2048 }
2049 }
2050
2051
2052 inline void
2053 PreconditionBase::vmult(MPI::Vector &dst, const MPI::Vector &src) const
2054 {
2055 Assert(dst.trilinos_partitioner().SameAs(
2056 preconditioner->OperatorRangeMap()),
2057 ExcNonMatchingMaps("dst"));
2058 Assert(src.trilinos_partitioner().SameAs(
2059 preconditioner->OperatorDomainMap()),
2060 ExcNonMatchingMaps("src"));
2061
2062 const int ierr = preconditioner->ApplyInverse(src.trilinos_vector(),
2063 dst.trilinos_vector());
2064 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2065 }
2066
2067 inline void
2068 PreconditionBase::Tvmult(MPI::Vector &dst, const MPI::Vector &src) const
2069 {
2070 Assert(dst.trilinos_partitioner().SameAs(
2071 preconditioner->OperatorRangeMap()),
2072 ExcNonMatchingMaps("dst"));
2073 Assert(src.trilinos_partitioner().SameAs(
2074 preconditioner->OperatorDomainMap()),
2075 ExcNonMatchingMaps("src"));
2076
2077 preconditioner->SetUseTranspose(true);
2078 const int ierr = preconditioner->ApplyInverse(src.trilinos_vector(),
2079 dst.trilinos_vector());
2080 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2081 preconditioner->SetUseTranspose(false);
2082 }
2083
2084
2085 // For the implementation of the <code>vmult</code> function with deal.II
2086 // data structures we note that invoking a call of the Trilinos
2087 // preconditioner requires us to use Epetra vectors as well. We do this by
2088 // providing a view, i.e., feed Trilinos with a pointer to the data, so we
2089 // avoid copying the content of the vectors during the iteration (this
2090 // function is only useful when used in serial anyway). In the declaration
2091 // of the right hand side, we need to cast the source vector (that is
2092 // <code>const</code> in all deal.II calls) to non-constant value, as this
2093 // is the way Trilinos wants to have them.
2094 inline void
2096 const ::Vector<double> &src) const
2097 {
2098 AssertDimension(dst.size(),
2099 preconditioner->OperatorDomainMap().NumMyElements());
2100 AssertDimension(src.size(),
2101 preconditioner->OperatorRangeMap().NumMyElements());
2102 Epetra_Vector tril_dst(View,
2103 preconditioner->OperatorDomainMap(),
2104 dst.begin());
2105 Epetra_Vector tril_src(View,
2106 preconditioner->OperatorRangeMap(),
2107 const_cast<double *>(src.begin()));
2108
2109 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2110 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2111 }
2112
2113
2114 inline void
2116 const ::Vector<double> &src) const
2117 {
2118 AssertDimension(dst.size(),
2119 preconditioner->OperatorDomainMap().NumMyElements());
2120 AssertDimension(src.size(),
2121 preconditioner->OperatorRangeMap().NumMyElements());
2122 Epetra_Vector tril_dst(View,
2123 preconditioner->OperatorDomainMap(),
2124 dst.begin());
2125 Epetra_Vector tril_src(View,
2126 preconditioner->OperatorRangeMap(),
2127 const_cast<double *>(src.begin()));
2128
2129 preconditioner->SetUseTranspose(true);
2130 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2131 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2132 preconditioner->SetUseTranspose(false);
2133 }
2134
2135
2136
2137 inline void
2141 {
2143 preconditioner->OperatorDomainMap().NumMyElements());
2145 preconditioner->OperatorRangeMap().NumMyElements());
2146 Epetra_Vector tril_dst(View,
2147 preconditioner->OperatorDomainMap(),
2148 dst.begin());
2149 Epetra_Vector tril_src(View,
2150 preconditioner->OperatorRangeMap(),
2151 const_cast<double *>(src.begin()));
2152
2153 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2154 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2155 }
2156
2157 inline void
2161 {
2163 preconditioner->OperatorDomainMap().NumMyElements());
2165 preconditioner->OperatorRangeMap().NumMyElements());
2166 Epetra_Vector tril_dst(View,
2167 preconditioner->OperatorDomainMap(),
2168 dst.begin());
2169 Epetra_Vector tril_src(View,
2170 preconditioner->OperatorRangeMap(),
2171 const_cast<double *>(src.begin()));
2172
2173 preconditioner->SetUseTranspose(true);
2174 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2175 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2176 preconditioner->SetUseTranspose(false);
2177 }
2178
2179# endif
2180
2181} // namespace TrilinosWrappers
2182
2183
2188
2189# endif // DEAL_II_WITH_TRILINOS
2190
2191#endif // trilinos_precondition_h
2192/*------------------------- trilinos_precondition.h -------------------------*/
Definition: vector.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583
void vmult(LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const override
void vmult(MPI::Vector &dst, const MPI::Vector &src) const override
std::shared_ptr< SparseMatrix > trilinos_matrix
virtual void Tvmult(::LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
virtual void Tvmult(MPI::Vector &dst, const MPI::Vector &src) const
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
~PreconditionBase() override=default
virtual ~PreconditionAMGMueLu() override=default
std::shared_ptr< SparseMatrix > trilinos_matrix
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
virtual void vmult(MPI::Vector &dst, const MPI::Vector &src) const
Epetra_Operator & trilinos_operator() const
static ::ExceptionBase & ExcNonMatchingMaps(std::string arg1)
std::shared_ptr< Epetra_Map > vector_distributor
Teuchos::RCP< Epetra_Operator > preconditioner
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
virtual void Tvmult(::Vector< double > &dst, const ::Vector< double > &src) 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 Tvmult(MPI::Vector &dst, const MPI::Vector &src) const override
virtual void vmult(::Vector< double > &dst, const ::Vector< double > &src) const
void Tvmult(LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const override
void set_operator_null_space(Teuchos::ParameterList &parameter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const
virtual void vmult(::LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) 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 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 initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
void set_parameters(Teuchos::ParameterList &parameter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const
size_type locally_owned_size() const
size_type size() const
iterator begin()
unsigned int global_dof_index
Definition: types.h:76