10#ifndef IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
11#define IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
13#include "Ifpack2_Factory.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"
33#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
34#include "Ifpack2_Details_Fic.hpp"
35#include "Ifpack2_Details_Fildl.hpp"
36#include "Ifpack2_Details_Filu.hpp"
39#ifdef HAVE_IFPACK2_AMESOS2
40#include "Ifpack2_Details_Amesos2Wrapper.hpp"
43#ifdef HAVE_IFPACK2_HYPRE
44#include "Ifpack2_Hypre.hpp"
50template <
class MatrixType>
51Teuchos::RCP<typename OneLevelFactory<MatrixType>::prec_type>
53 const Teuchos::RCP<const row_matrix_type>&
matrix)
const {
65 prec =
rcp(new ::Ifpack2::Chebyshev<row_matrix_type>(
matrix));
69#ifdef HAVE_IFPACK2_AMESOS2
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.");
91#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
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");
111 "The \"KRYLOV\" preconditioner option has "
112 "been deprecated and removed. If you want a Krylov solver, use the "
124 Teuchos::ParameterList
params;
125 params.set(
"relaxation: container",
"Dense");
138#ifdef HAVE_IFPACK2_AMESOS2
140 Teuchos::ParameterList
params;
141 params.set(
"relaxation: container",
"SparseAmesos2");
145 "Ifpack2::Details::OneLevelFactory: "
146 "\"SPARSE BLOCK RELAXATION\" requires building Trilinos with Amesos2 enabled.");
155 Teuchos::ParameterList
params;
156 params.set(
"relaxation: container",
"TriDi");
162 Teuchos::ParameterList
params;
163 params.set(
"relaxation: container",
"Banded");
169 else if (
precTypeUpper ==
"LOCAL SPARSE TRIANGULAR SOLVER" ||
179#ifdef HAVE_IFPACK2_HYPRE
186 true, std::invalid_argument,
187 "Ifpack2::Details::OneLevelFactory::create: "
188 "Invalid preconditioner type \""
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.");
200template <
class MatrixType>
201std::vector<std::string>
204 "CHEBYSHEV",
"DENSE",
"LAPACK",
205#ifdef HAVE_IFPACK2_AMESOS2
208 "DIAGONAL",
"ILUT",
"RELAXATION",
"RILUK",
"RBILUK",
"MDF",
209#ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
210 "FAST_IC",
"FAST_ILU",
"FAST_ILU_B",
"FAST_ILDL",
212 "BLOCK_RELAXATION",
"BLOCK RELAXATION",
"BLOCKRELAXATION",
"DENSE_BLOCK_RELAXATION",
"DENSE BLOCK RELAXATION",
"DENSEBLOCKRELAXATION",
214#ifdef HAVE_IFPACK2_AMESOS2
215 "SPARSE_BLOCK_RELAXATION",
"SPARSE BLOCK RELAXATION",
"SPARSEBLOCKRELAXATION",
217#ifdef HAVE_IFPACK2_HYPRE
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",
225 return supportedNames;
228template <
class MatrixType>
229bool OneLevelFactory<MatrixType>::isSupported(
const std::string& precType)
const {
231 std::string precTypeUpper = canonicalize(precType);
232 std::vector<std::string> supportedNames = getSupportedNames();
235 auto it = std::find(std::begin(supportedNames), std::end(supportedNames), precTypeUpper);
236 return it != std::end(supportedNames);
242#define IFPACK2_DETAILS_ONELEVELFACTORY_INSTANT(S, LO, GO, N) \
243 template class Ifpack2::Details::OneLevelFactory<Tpetra::RowMatrix<S, LO, GO, N> >;
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 > ¶ms)
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