43 #ifndef IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
44 #define IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
46 #include "Ifpack2_Chebyshev.hpp"
47 #include "Ifpack2_Details_DenseSolver.hpp"
48 #include "Ifpack2_Diagonal.hpp"
49 #include "Ifpack2_IdentitySolver.hpp"
50 #include "Ifpack2_ILUT.hpp"
51 #include "Ifpack2_Relaxation.hpp"
52 #include "Ifpack2_RILUK.hpp"
53 #include "Ifpack2_Experimental_RBILUK.hpp"
54 #include "Ifpack2_BlockRelaxation.hpp"
55 #include "Ifpack2_BandedContainer.hpp"
56 #include "Ifpack2_DenseContainer.hpp"
57 #include "Ifpack2_SparseContainer.hpp"
58 #include "Ifpack2_TriDiContainer.hpp"
59 #include "Ifpack2_LocalSparseTriangularSolver.hpp"
61 #ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
62 #include "Ifpack2_Details_Fic.hpp"
63 #include "Ifpack2_Details_Fildl.hpp"
64 #include "Ifpack2_Details_Filu.hpp"
65 #endif // HAVE_IFPACK2_SHYLU_NODEFASTILU
67 #ifdef HAVE_IFPACK2_AMESOS2
68 # include "Ifpack2_Details_Amesos2Wrapper.hpp"
69 #endif // HAVE_IFPACK2_AMESOS2
74 template<
class MatrixType>
85 std::string precTypeUpper (precType);
86 if (precTypeUpper.size () > 0) {
88 for (
size_t k = 0; k < precTypeUpper.size (); ++k) {
89 precTypeUpper[k] = std::toupper<char> (precTypeUpper[k], locale);
93 if (precTypeUpper ==
"CHEBYSHEV") {
96 prec =
rcp (new ::Ifpack2::Chebyshev<row_matrix_type> (matrix));
98 else if (precTypeUpper ==
"DENSE" || precTypeUpper ==
"LAPACK") {
101 else if (precTypeUpper ==
"AMESOS2") {
102 #ifdef HAVE_IFPACK2_AMESOS2
106 true, std::invalid_argument,
"Ifpack2::Details::OneLevelFactory: "
107 "You may not ask for the preconditioner \"AMESOS2\" unless "
108 "you have built Trilinos with the Amesos2 package enabled.");
109 #endif // HAVE_IFPACK2_AMESOS2
111 else if (precTypeUpper ==
"DIAGONAL") {
112 prec =
rcp (
new Diagonal<row_matrix_type> (matrix));
114 else if (precTypeUpper ==
"ILUT") {
117 else if (precTypeUpper ==
"RELAXATION") {
120 else if (precTypeUpper ==
"RILUK") {
123 else if (precTypeUpper ==
"RBILUK") {
126 else if (precTypeUpper ==
"FAST_IC" || precTypeUpper ==
"FAST_ILU" || precTypeUpper ==
"FAST_ILDL") {
127 #ifdef HAVE_IFPACK2_SHYLU_NODEFASTILU
129 if(precTypeUpper ==
"FAST_IC")
131 else if(precTypeUpper ==
"FAST_ILU")
133 else if(precTypeUpper ==
"FAST_ILDL")
138 throw std::invalid_argument(
"The Ifpack2 FastIC, FastILU and FastILDL preconditioners require the FastILU subpackage of ShyLU to be enabled\n"
139 "To enable FastILU, set the CMake option Trilinos_ENABLE_ShyLU_NodeFastILU=ON");
143 else if (precTypeUpper ==
"KRYLOV") {
145 (
true, std::invalid_argument,
"The \"KRYLOV\" preconditioner option has "
146 "been deprecated and removed. If you want a Krylov solver, use the "
149 else if (precTypeUpper ==
"BLOCK_RELAXATION" ||
150 precTypeUpper ==
"BLOCK RELAXATION" ||
151 precTypeUpper ==
"BLOCKRELAXATION" ||
152 precTypeUpper ==
"DENSE_BLOCK_RELAXATION" ||
153 precTypeUpper ==
"DENSE BLOCK RELAXATION" ||
154 precTypeUpper ==
"DENSEBLOCKRELAXATION" ) {
160 params.
set (
"relaxation: container",
"Dense");
163 else if (precTypeUpper ==
"SPARSE_BLOCK_RELAXATION" ||
164 precTypeUpper ==
"SPARSE BLOCK RELAXATION" ||
165 precTypeUpper ==
"SPARSEBLOCKRELAXATION" ) {
172 #ifdef HAVE_IFPACK2_AMESOS2
175 params.
set (
"relaxation: container",
"SparseAmesos2");
179 (
true, std::invalid_argument,
"Ifpack2::Details::OneLevelFactory: "
180 "\"SPARSE BLOCK RELAXATION\" requires building Trilinos with Amesos2 enabled.");
183 else if (precTypeUpper ==
"TRIDI_RELAXATION" ||
184 precTypeUpper ==
"TRIDI RELAXATION" ||
185 precTypeUpper ==
"TRIDIRELAXATION" ||
186 precTypeUpper ==
"TRIDIAGONAL_RELAXATION" ||
187 precTypeUpper ==
"TRIDIAGONAL RELAXATION" ||
188 precTypeUpper ==
"TRIDIAGONALRELAXATION") {
191 params.
set (
"relaxation: container",
"TriDi");
194 else if (precTypeUpper ==
"BANDED_RELAXATION" ||
195 precTypeUpper ==
"BANDED RELAXATION" ||
196 precTypeUpper ==
"BANDEDRELAXATION") {
199 params.
set (
"relaxation: container",
"Banded");
202 else if (precTypeUpper ==
"IDENTITY" || precTypeUpper ==
"IDENTITY_SOLVER") {
206 else if (precTypeUpper ==
"LOCAL SPARSE TRIANGULAR SOLVER" ||
207 precTypeUpper ==
"LOCAL_SPARSE_TRIANGULAR_SOLVER" ||
208 precTypeUpper ==
"LOCALSPARSETRIANGULARSOLVER" ||
209 precTypeUpper ==
"SPARSE TRIANGULAR SOLVER" ||
210 precTypeUpper ==
"SPARSE_TRIANGULAR_SOLVER" ||
211 precTypeUpper ==
"SPARSETRIANGULARSOLVER") {
216 true, std::invalid_argument,
"Ifpack2::Details::OneLevelFactory::create: "
217 "Invalid preconditioner type \"" << precType <<
"\".");
221 prec.is_null (), std::logic_error,
"Ifpack2::Details::OneLevelFactory::"
222 "create: Return value is null right before return. This should never "
223 "happen. Please report this bug to the Ifpack2 developers.");
230 #define IFPACK2_DETAILS_ONELEVELFACTORY_INSTANT(S,LO,GO,N) \
231 template class Ifpack2::Details::OneLevelFactory< Tpetra::RowMatrix<S, LO, GO, N> >;
233 #endif // IFPACK2_DETAILS_ONELEVELFACTORY_DEF_HPP
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:76
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:243
"Preconditioner" that uses LAPACK's dense LU.
Definition: Ifpack2_Details_DenseSolver_decl.hpp:73
The Ifpack2 wrapper to the ILDL preconditioner of ShyLU FastILU.
Definition: Ifpack2_Details_Fildl_decl.hpp:62
The Ifpack2 wrapper to the ILU preconditioner of ShyLU FastILU.
Definition: Ifpack2_Details_Filu_decl.hpp:62
ILUT (incomplete LU factorization with threshold) of a Tpetra sparse matrix.
Definition: Ifpack2_ILUT_decl.hpp:91
Block relaxation preconditioners (or smoothers) for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse ma...
Definition: Ifpack2_BlockRelaxation_decl.hpp:83
"Identity" preconditioner.
Definition: Ifpack2_IdentitySolver_decl.hpp:59
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Wrapper class for direct solvers in Amesos2.
Definition: Ifpack2_Details_Amesos2Wrapper_decl.hpp:102
"Preconditioner" that solves local sparse triangular systems.
Definition: Ifpack2_LocalSparseTriangularSolver_decl.hpp:77
ParameterList & setParameters(const ParameterList &source)
The Ifpack2 wrapper to the incomplete Chebyshev preconditioner of ShyLU FastILU.
Definition: Ifpack2_Details_Fic_decl.hpp:62
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:223
ILU(k) factorization of a given Tpetra::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:128