deal.II version 9.7.0
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2008 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_trilinos_precondition_h
16#define dealii_trilinos_precondition_h
17
18
19#include <deal.II/base/config.h>
20
21#ifdef DEAL_II_WITH_TRILINOS
22
24
27
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>
36
37# include <memory>
38
39// forward declarations
40# ifndef DOXYGEN
41class Ifpack_Preconditioner;
42class Ifpack_Chebyshev;
43namespace ML_Epetra
44{
45 class MultiLevelPreconditioner;
46}
47# endif
48
50
51// forward declarations
52# ifndef DOXYGEN
53template <typename number>
54class SparseMatrix;
55template <typename number>
56class Vector;
57class SparsityPattern;
58# endif
59
64
65namespace TrilinosWrappers
66{
67 // forward declarations
68 class SparseMatrix;
70 class SolverBase;
71
79 {
80 public:
85
91 {};
92
99
104 void
105 clear();
106
111 get_mpi_communicator() const;
112
122 void
124
128 virtual void
129 vmult(MPI::Vector &dst, const MPI::Vector &src) const;
130
134 virtual void
135 Tvmult(MPI::Vector &dst, const MPI::Vector &src) const;
136
141 virtual void
142 vmult(::Vector<double> &dst, const ::Vector<double> &src) const;
143
148 virtual void
150 const ::Vector<double> &src) const;
151
156 virtual void
158 const ::LinearAlgebra::distributed::Vector<double> &src) const;
159
164 virtual void
166 const ::LinearAlgebra::distributed::Vector<double> &src) const;
167
178 trilinos_operator() const;
180
185
192
200
202
211 std::string,
212 << "The sparse matrix the preconditioner is based on "
213 << "uses a map that is not compatible to the one in vector "
214 << arg1 << ". Check preconditioner and matrix setup.");
216
217 friend class SolverBase;
218
219 protected:
224 Teuchos::RCP<Epetra_Operator> preconditioner;
225
230 Epetra_MpiComm communicator;
231 };
232
233
250 {
251 public:
264 {
269 AdditionalData(const double omega = 1,
270 const double min_diagonal = 0,
271 const unsigned int n_sweeps = 1);
272
276 double omega;
277
286
291 unsigned int n_sweeps;
292 };
293
298 void
299 initialize(const SparseMatrix &matrix,
300 const AdditionalData &additional_data = AdditionalData());
301 };
302
303
304
331 {
332 public:
347 {
354 AdditionalData(const double omega = 1,
355 const double min_diagonal = 0,
356 const unsigned int overlap = 0,
357 const unsigned int n_sweeps = 1);
358
363 double omega;
364
373
378 unsigned int overlap;
379
384 unsigned int n_sweeps;
385 };
386
392 void
393 initialize(const SparseMatrix &matrix,
394 const AdditionalData &additional_data = AdditionalData());
395 };
396
397
398
425 {
426 public:
441 {
448 AdditionalData(const double omega = 1,
449 const double min_diagonal = 0,
450 const unsigned int overlap = 0,
451 const unsigned int n_sweeps = 1);
452
457 double omega;
458
467
472 unsigned int overlap;
473
478 unsigned int n_sweeps;
479 };
480
486 void
487 initialize(const SparseMatrix &matrix,
488 const AdditionalData &additional_data = AdditionalData());
489 };
490
491
492
510 {
511 public:
528 {
534 AdditionalData(const unsigned int block_size = 1,
535 const std::string &block_creation_type = "linear",
536 const double omega = 1,
537 const double min_diagonal = 0,
538 const unsigned int n_sweeps = 1);
539
543 unsigned int block_size;
544
553
558 double omega;
559
568
573 unsigned int n_sweeps;
574 };
575
580 void
581 initialize(const SparseMatrix &matrix,
582 const AdditionalData &additional_data = AdditionalData());
583 };
584
585
586
609 {
610 public:
629 {
637 AdditionalData(const unsigned int block_size = 1,
638 const std::string &block_creation_type = "linear",
639 const double omega = 1,
640 const double min_diagonal = 0,
641 const unsigned int overlap = 0,
642 const unsigned int n_sweeps = 1);
643
647 unsigned int block_size;
648
657
662 double omega;
663
672
677 unsigned int overlap;
678
683 unsigned int n_sweeps;
684 };
685
691 void
692 initialize(const SparseMatrix &matrix,
693 const AdditionalData &additional_data = AdditionalData());
694 };
695
696
697
720 {
721 public:
740 {
748 AdditionalData(const unsigned int block_size = 1,
749 const std::string &block_creation_type = "linear",
750 const double omega = 1,
751 const double min_diagonal = 0,
752 const unsigned int overlap = 0,
753 const unsigned int n_sweeps = 1);
754
758 unsigned int block_size;
759
768
773 double omega;
774
783
788 unsigned int overlap;
789
794 unsigned int n_sweeps;
795 };
796
802 void
803 initialize(const SparseMatrix &matrix,
804 const AdditionalData &additional_data = AdditionalData());
805 };
806
807
808
844 {
845 public:
863 {
873 AdditionalData(const unsigned int ic_fill = 0,
874 const double ic_atol = 0.,
875 const double ic_rtol = 1.,
876 const unsigned int overlap = 0);
877
886 unsigned int ic_fill;
887
893 double ic_atol;
894
899 double ic_rtol;
900
905 unsigned int overlap;
906 };
907
912 void
913 initialize(const SparseMatrix &matrix,
914 const AdditionalData &additional_data = AdditionalData());
915 };
916
917
918
948 {
949 public:
988 {
992 AdditionalData(const unsigned int ilu_fill = 0,
993 const double ilu_atol = 0.,
994 const double ilu_rtol = 1.,
995 const unsigned int overlap = 0);
996
1000 unsigned int ilu_fill;
1001
1006 double ilu_atol;
1007
1012 double ilu_rtol;
1013
1017 unsigned int overlap;
1018 };
1019
1024 void
1025 initialize(const SparseMatrix &matrix,
1026 const AdditionalData &additional_data = AdditionalData());
1027 };
1028
1029
1030
1066 {
1067 public:
1086 {
1097 AdditionalData(const double ilut_drop = 0.,
1098 const unsigned int ilut_fill = 0,
1099 const double ilut_atol = 0.,
1100 const double ilut_rtol = 1.,
1101 const unsigned int overlap = 0);
1102
1108
1117 unsigned int ilut_fill;
1118
1125
1131
1136 unsigned int overlap;
1137 };
1138
1143 void
1144 initialize(const SparseMatrix &matrix,
1145 const AdditionalData &additional_data = AdditionalData());
1146 };
1147
1148
1149
1168 {
1169 public:
1175 {
1179 AdditionalData(const unsigned int overlap = 0);
1180
1181
1186 unsigned int overlap;
1187 };
1188
1193 void
1194 initialize(const SparseMatrix &matrix,
1195 const AdditionalData &additional_data = AdditionalData());
1196 };
1197
1198
1199
1209 {
1210 public:
1216 {
1220 AdditionalData(const unsigned int degree = 1,
1221 const double max_eigenvalue = 10.,
1222 const double eigenvalue_ratio = 30.,
1223 const double min_eigenvalue = 1.,
1224 const double min_diagonal = 1e-12,
1225 const bool nonzero_starting = false);
1226
1232 unsigned int degree;
1233
1239
1244
1250
1256
1268 };
1269
1274 void
1275 initialize(const SparseMatrix &matrix,
1276 const AdditionalData &additional_data = AdditionalData());
1277 };
1278
1279
1280
1323 {
1324 public:
1332 {
1356 AdditionalData(const bool elliptic = true,
1357 const bool higher_order_elements = false,
1358 const unsigned int n_cycles = 1,
1359 const bool w_cycle = false,
1360 const double aggregation_threshold = 1e-4,
1361 const std::vector<std::vector<bool>> &constant_modes =
1362 std::vector<std::vector<bool>>(0),
1363 const unsigned int smoother_sweeps = 2,
1364 const unsigned int smoother_overlap = 0,
1365 const bool output_details = false,
1366 const char *smoother_type = "Chebyshev",
1367 const char *coarse_type = "Amesos-KLU");
1368
1399 void
1401 Teuchos::ParameterList &parameter_list,
1402 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1403 const Epetra_RowMatrix &matrix) const;
1404
1412 void
1414 Teuchos::ParameterList &parameter_list,
1415 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1416 const SparseMatrix &matrix) const;
1417
1423 void
1425 Teuchos::ParameterList &parameter_list,
1426 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1427 const Epetra_RowMatrix &matrix) const;
1428
1434 void
1436 Teuchos::ParameterList &parameter_list,
1437 std::unique_ptr<Epetra_MultiVector> &distributed_constant_modes,
1438 const SparseMatrix &matrix) const;
1439
1448
1454
1459 unsigned int n_cycles;
1460
1466
1477
1502 std::vector<std::vector<bool>> constant_modes;
1503
1511 std::vector<std::vector<double>> constant_modes_values;
1512
1523 unsigned int smoother_sweeps;
1524
1529 unsigned int smoother_overlap;
1530
1537
1572 const char *smoother_type;
1573
1578 const char *coarse_type;
1579 };
1580
1584 ~PreconditionAMG() override;
1585
1586
1592 void
1593 initialize(const SparseMatrix &matrix,
1594 const AdditionalData &additional_data = AdditionalData());
1595
1614 void
1615 initialize(const Epetra_RowMatrix &matrix,
1616 const AdditionalData &additional_data = AdditionalData());
1617
1632 void
1633 initialize(const SparseMatrix &matrix,
1634 const Teuchos::ParameterList &ml_parameters);
1635
1643 void
1644 initialize(const Epetra_RowMatrix &matrix,
1645 const Teuchos::ParameterList &ml_parameters);
1646
1653 template <typename number>
1654 void
1655 initialize(const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1656 const AdditionalData &additional_data = AdditionalData(),
1657 const double drop_tolerance = 1e-13,
1658 const ::SparsityPattern *use_this_sparsity = nullptr);
1659
1672 void
1673 reinit();
1674
1679 void
1680 clear();
1681
1685 size_type
1686 memory_consumption() const;
1687
1688 private:
1692 std::shared_ptr<SparseMatrix> trilinos_matrix;
1693 };
1694
1695
1696
1697# if defined(DOXYGEN) || defined(DEAL_II_TRILINOS_WITH_MUELU)
1725 {
1726 public:
1734 {
1739 AdditionalData(const bool elliptic = true,
1740 const unsigned int n_cycles = 1,
1741 const bool w_cycle = false,
1742 const double aggregation_threshold = 1e-4,
1743 const std::vector<std::vector<bool>> &constant_modes =
1744 std::vector<std::vector<bool>>(0),
1745 const unsigned int smoother_sweeps = 2,
1746 const unsigned int smoother_overlap = 0,
1747 const bool output_details = false,
1748 const char *smoother_type = "Chebyshev",
1749 const char *coarse_type = "Amesos-KLU");
1750
1759
1764 unsigned int n_cycles;
1765
1771
1782
1789 std::vector<std::vector<bool>> constant_modes;
1790
1801 unsigned int smoother_sweeps;
1802
1807 unsigned int smoother_overlap;
1808
1815
1850 const char *smoother_type;
1851
1856 const char *coarse_type;
1857 };
1858
1863
1867 virtual ~PreconditionAMGMueLu() override = default;
1868
1874 void
1875 initialize(const SparseMatrix &matrix,
1876 const AdditionalData &additional_data = AdditionalData());
1877
1884 void
1885 initialize(const Epetra_CrsMatrix &matrix,
1886 const AdditionalData &additional_data = AdditionalData());
1887
1901 void
1902 initialize(const SparseMatrix &matrix,
1903 Teuchos::ParameterList &muelu_parameters);
1904
1910 void
1911 initialize(const Epetra_CrsMatrix &matrix,
1912 Teuchos::ParameterList &muelu_parameters);
1913
1920 template <typename number>
1921 void
1922 initialize(const ::SparseMatrix<number> &deal_ii_sparse_matrix,
1923 const AdditionalData &additional_data = AdditionalData(),
1924 const double drop_tolerance = 1e-13,
1925 const ::SparsityPattern *use_this_sparsity = nullptr);
1926
1931 void
1932 clear();
1933
1937 size_type
1938 memory_consumption() const;
1939
1940 private:
1944 std::shared_ptr<SparseMatrix> trilinos_matrix;
1945 };
1946# endif
1947
1948
1949
1957 {
1958 public:
1964 {};
1965
1972 void
1973 initialize(const SparseMatrix &matrix,
1974 const AdditionalData &additional_data = AdditionalData());
1975
1979 void
1980 vmult(MPI::Vector &dst, const MPI::Vector &src) const override;
1981
1985 void
1986 Tvmult(MPI::Vector &dst, const MPI::Vector &src) const override;
1987
1992 void
1994 const ::Vector<double> &src) const override;
1995
2000 void
2002 const ::Vector<double> &src) const override;
2003
2008 void
2010 const ::LinearAlgebra::distributed::Vector<double> &src)
2011 const override;
2012
2018 void
2020 const ::LinearAlgebra::distributed::Vector<double> &src)
2021 const override;
2022 };
2023
2024
2025
2026 // ----------------------- inline and template functions --------------------
2027
2028
2029# ifndef DOXYGEN
2030
2031
2032 inline void
2034 {
2035 // This only flips a flag that tells
2036 // Trilinos that any vmult operation
2037 // should be done with the
2038 // transpose. However, the matrix
2039 // structure is not reset.
2040 int ierr;
2041
2042 if (!preconditioner->UseTranspose())
2043 {
2044 ierr = preconditioner->SetUseTranspose(true);
2045 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2046 }
2047 else
2048 {
2049 ierr = preconditioner->SetUseTranspose(false);
2050 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2051 }
2052 }
2053
2054
2055 inline void
2056 PreconditionBase::vmult(MPI::Vector &dst, const MPI::Vector &src) const
2057 {
2058 Assert(dst.trilinos_partitioner().SameAs(
2059 preconditioner->OperatorRangeMap()),
2060 ExcNonMatchingMaps("dst"));
2061 Assert(src.trilinos_partitioner().SameAs(
2062 preconditioner->OperatorDomainMap()),
2063 ExcNonMatchingMaps("src"));
2064
2065 const int ierr = preconditioner->ApplyInverse(src.trilinos_vector(),
2066 dst.trilinos_vector());
2067 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2068 }
2069
2070 inline void
2071 PreconditionBase::Tvmult(MPI::Vector &dst, const MPI::Vector &src) const
2072 {
2073 Assert(dst.trilinos_partitioner().SameAs(
2074 preconditioner->OperatorRangeMap()),
2075 ExcNonMatchingMaps("dst"));
2076 Assert(src.trilinos_partitioner().SameAs(
2077 preconditioner->OperatorDomainMap()),
2078 ExcNonMatchingMaps("src"));
2079
2080 preconditioner->SetUseTranspose(true);
2081 const int ierr = preconditioner->ApplyInverse(src.trilinos_vector(),
2082 dst.trilinos_vector());
2083 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2084 preconditioner->SetUseTranspose(false);
2085 }
2086
2087
2088 // For the implementation of the <code>vmult</code> function with deal.II
2089 // data structures we note that invoking a call of the Trilinos
2090 // preconditioner requires us to use Epetra vectors as well. We do this by
2091 // providing a view, i.e., feed Trilinos with a pointer to the data, so we
2092 // avoid copying the content of the vectors during the iteration (this
2093 // function is only useful when used in serial anyway). In the declaration
2094 // of the right hand side, we need to cast the source vector (that is
2095 // <code>const</code> in all deal.II calls) to non-constant value, as this
2096 // is the way Trilinos wants to have them.
2097 inline void
2098 PreconditionBase::vmult(::Vector<double> &dst,
2099 const ::Vector<double> &src) const
2100 {
2101 AssertDimension(dst.size(),
2102 preconditioner->OperatorDomainMap().NumMyElements());
2103 AssertDimension(src.size(),
2104 preconditioner->OperatorRangeMap().NumMyElements());
2105 Epetra_Vector tril_dst(View,
2106 preconditioner->OperatorDomainMap(),
2107 dst.begin());
2108 Epetra_Vector tril_src(View,
2109 preconditioner->OperatorRangeMap(),
2110 const_cast<double *>(src.begin()));
2111
2112 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2113 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2114 }
2115
2116
2117 inline void
2118 PreconditionBase::Tvmult(::Vector<double> &dst,
2119 const ::Vector<double> &src) const
2120 {
2121 AssertDimension(dst.size(),
2122 preconditioner->OperatorDomainMap().NumMyElements());
2123 AssertDimension(src.size(),
2124 preconditioner->OperatorRangeMap().NumMyElements());
2125 Epetra_Vector tril_dst(View,
2126 preconditioner->OperatorDomainMap(),
2127 dst.begin());
2128 Epetra_Vector tril_src(View,
2129 preconditioner->OperatorRangeMap(),
2130 const_cast<double *>(src.begin()));
2131
2132 preconditioner->SetUseTranspose(true);
2133 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2134 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2135 preconditioner->SetUseTranspose(false);
2136 }
2137
2138
2139
2140 inline void
2142 LinearAlgebra::distributed::Vector<double> &dst,
2143 const LinearAlgebra::distributed::Vector<double> &src) const
2144 {
2146 preconditioner->OperatorDomainMap().NumMyElements());
2148 preconditioner->OperatorRangeMap().NumMyElements());
2149 Epetra_Vector tril_dst(View,
2150 preconditioner->OperatorDomainMap(),
2151 dst.begin());
2152 Epetra_Vector tril_src(View,
2153 preconditioner->OperatorRangeMap(),
2154 const_cast<double *>(src.begin()));
2155
2156 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2157 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2158 }
2159
2160 inline void
2162 LinearAlgebra::distributed::Vector<double> &dst,
2163 const LinearAlgebra::distributed::Vector<double> &src) const
2164 {
2166 preconditioner->OperatorDomainMap().NumMyElements());
2168 preconditioner->OperatorRangeMap().NumMyElements());
2169 Epetra_Vector tril_dst(View,
2170 preconditioner->OperatorDomainMap(),
2171 dst.begin());
2172 Epetra_Vector tril_src(View,
2173 preconditioner->OperatorRangeMap(),
2174 const_cast<double *>(src.begin()));
2175
2176 preconditioner->SetUseTranspose(true);
2177 const int ierr = preconditioner->ApplyInverse(tril_src, tril_dst);
2178 AssertThrow(ierr == 0, ExcTrilinosError(ierr));
2179 preconditioner->SetUseTranspose(false);
2180 }
2181
2182# endif
2183
2184} // namespace TrilinosWrappers
2185
2186
2188
2189
2191
2192#else
2193
2194// Make sure the scripts that create the C++20 module input files have
2195// something to latch on if the preprocessor #ifdef above would
2196// otherwise lead to an empty content of the file.
2199
2200#endif // DEAL_II_WITH_TRILINOS
2201
2202#endif
size_type locally_owned_size() const
std::shared_ptr< SparseMatrix > trilinos_matrix
virtual ~PreconditionAMGMueLu() override=default
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
std::shared_ptr< SparseMatrix > trilinos_matrix
void initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
virtual void Tvmult(::LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const
virtual void Tvmult(MPI::Vector &dst, const MPI::Vector &src) const
virtual void vmult(MPI::Vector &dst, const MPI::Vector &src) const
Epetra_Operator & trilinos_operator() const
Teuchos::RCP< Epetra_Operator > preconditioner
virtual void Tvmult(::Vector< double > &dst, const ::Vector< double > &src) const
virtual void vmult(::Vector< double > &dst, const ::Vector< double > &src) 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 initialize(const SparseMatrix &matrix, const AdditionalData &additional_data=AdditionalData())
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
void Tvmult(MPI::Vector &dst, const MPI::Vector &src) const override
void Tvmult(LinearAlgebra::distributed::Vector< double > &dst, const ::LinearAlgebra::distributed::Vector< double > &src) const override
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())
virtual size_type size() const override
iterator begin()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
static ::ExceptionBase & ExcNonMatchingMaps(std::string arg1)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcTrilinosError(int arg1)
#define AssertThrow(cond, exc)
unsigned int global_dof_index
Definition types.h:94
AdditionalData(const bool elliptic=true, const unsigned int n_cycles=1, const bool w_cycle=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 bool elliptic=true, const bool higher_order_elements=false, const unsigned int n_cycles=1, const bool w_cycle=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::vector< std::vector< double > > constant_modes_values
void set_operator_null_space(Teuchos::ParameterList &parameter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const
void set_parameters(Teuchos::ParameterList &parameter_list, std::unique_ptr< Epetra_MultiVector > &distributed_constant_modes, const Epetra_RowMatrix &matrix) const
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)
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)
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)
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 unsigned int ic_fill=0, const double ic_atol=0., const double ic_rtol=1., const unsigned int overlap=0)
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 unsigned int ilu_fill=0, const double ilu_atol=0., const double ilu_rtol=1., const unsigned int overlap=0)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int n_sweeps=1)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)
AdditionalData(const double omega=1, const double min_diagonal=0, const unsigned int overlap=0, const unsigned int n_sweeps=1)