Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_OneLevelFactory_def.hpp
1// @HEADER
2// *****************************************************************************
3// Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4//
5// Copyright 2009 NTESS and the Ifpack2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
11#define IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
12
13#include "Ifpack2_Factory.hpp"
14#include "Ifpack2_Utilities.hpp"
15#include "Ifpack2_Chebyshev.hpp"
16#include "Ifpack2_Details_DenseSolver.hpp"
17#include "Ifpack2_Diagonal.hpp"
18#include "Ifpack2_IdentitySolver.hpp"
19#include "Ifpack2_ILUT.hpp"
20#include "Ifpack2_MDF.hpp"
21#include "Ifpack2_Relaxation.hpp"
22#include "Ifpack2_RILUK.hpp"
23#include "Ifpack2_Experimental_RBILUK.hpp"
24#include "Ifpack2_BlockRelaxation.hpp"
25#include "Ifpack2_BandedContainer.hpp"
26#include "Ifpack2_DenseContainer.hpp"
27#include "Ifpack2_DatabaseSchwarz.hpp"
28#include "Ifpack2_SparseContainer.hpp"
29#include "Ifpack2_TriDiContainer.hpp"
30#include "Ifpack2_LocalSparseTriangularSolver.hpp"
31#include "Ifpack2_Hiptmair.hpp"
32
33#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
34#include "Ifpack2_Details_Fic.hpp"
35#include "Ifpack2_Details_Fildl.hpp"
36#include "Ifpack2_Details_Filu.hpp"
37#endif // HAVE_IFPACK2_SHYLU_NODEFASTILU
38
39#ifdef HAVE_IFPACK2_AMESOS2
40#include "Ifpack2_Details_Amesos2Wrapper.hpp"
41#endif // HAVE_IFPACK2_AMESOS2
42
43#ifdef HAVE_IFPACK2_HYPRE
44#include "Ifpack2_Hypre.hpp"
45#endif // HAVE_IFPACK2_HYPRE
46
47namespace Ifpack2 {
48namespace Details {
49
50template <class MatrixType>
51Teuchos::RCP<typename OneLevelFactory<MatrixType>::prec_type>
53 const Teuchos::RCP<const row_matrix_type>& matrix) const {
54 using Teuchos::RCP;
55 using Teuchos::rcp;
56
58
59 // precTypeUpper is the upper-case version of precType.
60 std::string precTypeUpper = canonicalize(precType);
61
62 if (precTypeUpper == "CHEBYSHEV") {
63 // We have to distinguish Ifpack2::Chebyshev from its
64 // implementation class Ifpack2::Details::Chebyshev.
65 prec = rcp(new ::Ifpack2::Chebyshev<row_matrix_type>(matrix));
66 } else if (precTypeUpper == "DENSE" || precTypeUpper == "LAPACK") {
68 } else if (precTypeUpper == "AMESOS2") {
69#ifdef HAVE_IFPACK2_AMESOS2
71#else
73 true, std::invalid_argument,
74 "Ifpack2::Details::OneLevelFactory: "
75 "You may not ask for the preconditioner \"AMESOS2\" unless "
76 "you have built Trilinos with the Amesos2 package enabled.");
77#endif // HAVE_IFPACK2_AMESOS2
78 } else if (precTypeUpper == "DIAGONAL") {
80 } else if (precTypeUpper == "ILUT") {
82 } else if (precTypeUpper == "RELAXATION") {
84 } else if (precTypeUpper == "RILUK") {
86 } else if (precTypeUpper == "MDF") {
88 } else if (precTypeUpper == "RBILUK") {
90 } else if (precTypeUpper == "FAST_IC" || precTypeUpper == "FAST_ILU" || precTypeUpper == "FAST_ILU_B" || precTypeUpper == "FAST_ILDL") {
91#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
92 {
93 if (precTypeUpper == "FAST_IC")
95 else if (precTypeUpper == "FAST_ILU")
97 else if (precTypeUpper == "FAST_ILU_B")
99 else if (precTypeUpper == "FAST_ILDL")
101 }
102#else
103 {
104 throw std::invalid_argument(
105 "The Ifpack2 FastIC, FastILU and FastILDL preconditioners require the FastILU subpackage of ShyLU to be enabled\n"
106 "To enable FastILU, set the CMake option Trilinos_ENABLE_ShyLU_NodeFastILU=ON");
107 }
108#endif
109 } else if (precTypeUpper == "KRYLOV") {
110 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
111 "The \"KRYLOV\" preconditioner option has "
112 "been deprecated and removed. If you want a Krylov solver, use the "
113 "Belos package.");
114 } else if (precTypeUpper == "BLOCK_RELAXATION" ||
115 precTypeUpper == "BLOCK RELAXATION" ||
116 precTypeUpper == "BLOCKRELAXATION" ||
117 precTypeUpper == "DENSE_BLOCK_RELAXATION" ||
118 precTypeUpper == "DENSE BLOCK RELAXATION" ||
119 precTypeUpper == "DENSEBLOCKRELAXATION") {
120 // NOTE (mfh 12 Aug 2016) Choice of "container type" is now a
121 // run-time parameter. The "ContainerType" template parameter is
122 // now always Container<row_matrix_type>.
124 Teuchos::ParameterList params;
125 params.set("relaxation: container", "Dense");
127 } else if (precTypeUpper == "DATABASE SCHWARZ") {
129 } else if (precTypeUpper == "SPARSE_BLOCK_RELAXATION" ||
130 precTypeUpper == "SPARSE BLOCK RELAXATION" ||
131 precTypeUpper == "SPARSEBLOCKRELAXATION") {
132 // FIXME (mfh 22 May 2014) We would prefer to have the choice of
133 // dense or sparse blocks (the "container type") be a run-time
134 // decision. This will require refactoring BlockRelaxation so
135 // that the "container type" is not a template parameter. For
136 // now, we default to use dense blocks.
137 // typedef SparseContainer<row_matrix_type, ILUT<row_matrix_type>> container_type;
138#ifdef HAVE_IFPACK2_AMESOS2
140 Teuchos::ParameterList params;
141 params.set("relaxation: container", "SparseAmesos2");
143#else
144 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument,
145 "Ifpack2::Details::OneLevelFactory: "
146 "\"SPARSE BLOCK RELAXATION\" requires building Trilinos with Amesos2 enabled.");
147#endif
148 } else if (precTypeUpper == "TRIDI_RELAXATION" ||
149 precTypeUpper == "TRIDI RELAXATION" ||
150 precTypeUpper == "TRIDIRELAXATION" ||
151 precTypeUpper == "TRIDIAGONAL_RELAXATION" ||
152 precTypeUpper == "TRIDIAGONAL RELAXATION" ||
153 precTypeUpper == "TRIDIAGONALRELAXATION") {
155 Teuchos::ParameterList params;
156 params.set("relaxation: container", "TriDi");
158 } else if (precTypeUpper == "BANDED_RELAXATION" ||
159 precTypeUpper == "BANDED RELAXATION" ||
160 precTypeUpper == "BANDEDRELAXATION") {
162 Teuchos::ParameterList params;
163 params.set("relaxation: container", "Banded");
165 } else if (precTypeUpper == "IDENTITY" || precTypeUpper == "IDENTITY_SOLVER") {
167 }
168
169 else if (precTypeUpper == "LOCAL SPARSE TRIANGULAR SOLVER" ||
170 precTypeUpper == "LOCAL_SPARSE_TRIANGULAR_SOLVER" ||
171 precTypeUpper == "LOCALSPARSETRIANGULARSOLVER" ||
172 precTypeUpper == "SPARSE TRIANGULAR SOLVER" ||
173 precTypeUpper == "SPARSE_TRIANGULAR_SOLVER" ||
174 precTypeUpper == "SPARSETRIANGULARSOLVER") {
176 } else if (precTypeUpper == "HIPTMAIR") {
178 }
179#ifdef HAVE_IFPACK2_HYPRE
180 else if (precTypeUpper == "HYPRE") {
182 }
183#endif
184 else {
186 true, std::invalid_argument,
187 "Ifpack2::Details::OneLevelFactory::create: "
188 "Invalid preconditioner type \""
189 << precType << "\".");
190 }
191
193 prec.is_null(), std::logic_error,
194 "Ifpack2::Details::OneLevelFactory::"
195 "create: Return value is null right before return. This should never "
196 "happen. Please report this bug to the Ifpack2 developers.");
197 return prec;
198}
199
200template <class MatrixType>
201std::vector<std::string>
203 std::vector<std::string> supportedNames = {
204 "CHEBYSHEV", "DENSE", "LAPACK",
205#ifdef HAVE_IFPACK2_AMESOS2
206 "AMESOS2",
207#endif
208 "DIAGONAL", "ILUT", "RELAXATION", "RILUK", "RBILUK", "MDF",
209#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
210 "FAST_IC", "FAST_ILU", "FAST_ILU_B", "FAST_ILDL",
211#endif
212 "BLOCK_RELAXATION", "BLOCK RELAXATION", "BLOCKRELAXATION", "DENSE_BLOCK_RELAXATION", "DENSE BLOCK RELAXATION", "DENSEBLOCKRELAXATION",
213 "DATABASE SCHWARZ",
214#ifdef HAVE_IFPACK2_AMESOS2
215 "SPARSE_BLOCK_RELAXATION", "SPARSE BLOCK RELAXATION", "SPARSEBLOCKRELAXATION",
216#endif
217#ifdef HAVE_IFPACK2_HYPRE
218 "HYPRE",
219#endif
220 "TRIDI_RELAXATION", "TRIDI RELAXATION", "TRIDIRELAXATION", "TRIDIAGONAL_RELAXATION", "TRIDIAGONAL RELAXATION", "TRIDIAGONALRELAXATION",
221 "BANDED_RELAXATION", "BANDED RELAXATION", "BANDEDRELAXATION",
222 "IDENTITY", "IDENTITY_SOLVER",
223 "LOCAL SPARSE TRIANGULAR SOLVER", "LOCAL_SPARSE_TRIANGULAR_SOLVER", "LOCALSPARSETRIANGULARSOLVER", "SPARSE TRIANGULAR SOLVER", "SPARSE_TRIANGULAR_SOLVER", "SPARSETRIANGULARSOLVER",
224 "HIPTMAIR"};
225 return supportedNames;
226}
227
228template <class MatrixType>
229bool OneLevelFactory<MatrixType>::isSupported(const std::string& precType) const {
230 // precTypeUpper is the upper-case version of precType.
231 std::string precTypeUpper = canonicalize(precType);
232 std::vector<std::string> supportedNames = getSupportedNames();
233 // const size_t numSupportedNames = supportedNames.size();
234 // const auto end = supportedNames + numSupportedNames;
235 auto it = std::find(std::begin(supportedNames), std::end(supportedNames), precTypeUpper);
236 return it != std::end(supportedNames);
237}
238
239} // namespace Details
240} // namespace Ifpack2
241
242#define IFPACK2_DETAILS_ONELEVELFACTORY_INSTANT(S, LO, GO, N) \
243 template class Ifpack2::Details::OneLevelFactory<Tpetra::RowMatrix<S, LO, GO, N> >;
244
245#endif // IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
File for utility functions.
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params)
Set the solver's parameters.
Definition Ifpack2_Details_LinearSolver_def.hpp:135
"Factory" for creating single-level preconditioners.
Definition Ifpack2_Details_OneLevelFactory_decl.hpp:92
Teuchos::RCP< prec_type > create(const std::string &precType, const Teuchos::RCP< const row_matrix_type > &matrix) const
Create an instance of Preconditioner given the string name of the preconditioner type.
Definition Ifpack2_Details_OneLevelFactory_def.hpp:52
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition Ifpack2_Experimental_RBILUK_decl.hpp:96
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40