46 #ifndef IFPACK2_CRSRILUK_DECL_HPP 
   47 #define IFPACK2_CRSRILUK_DECL_HPP 
   51 #include "Tpetra_CrsMatrix_decl.hpp" 
   54 #include "Ifpack2_LocalSparseTriangularSolver_decl.hpp" 
   56 #include <type_traits> 
  242 template<
class MatrixType>
 
  245                                            typename MatrixType::local_ordinal_type,
 
  246                                            typename MatrixType::global_ordinal_type,
 
  247                                            typename MatrixType::node_type>,
 
  249                                                                        typename MatrixType::local_ordinal_type,
 
  250                                                                        typename MatrixType::global_ordinal_type,
 
  251                                                                        typename MatrixType::node_type> >
 
  276   static_assert(std::is_same<MatrixType, row_matrix_type>::value, 
"Ifpack2::RILUK: The template parameter MatrixType must be a Tpetra::RowMatrix specialization.  Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.");
 
  284   template <
class NewMatrixType> 
friend class RILUK;
 
  311   template <
typename NewMatrixType>
 
  342     return isInitialized_;
 
  351     return numInitialize_;
 
  364     return initializeTime_;
 
  457   apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
 
  458          Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
 
  487   multiply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
 
  488             Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
 
  512       true, std::logic_error, 
"Ifpack2::RILUK::SetOverlapMode: " 
  513       "RILUK no longer implements overlap on its own.  " 
  514       "Use RILUK with AdditiveSchwarz if you want overlap.");
 
  519     return getL ().getGlobalNumEntries () + 
getU ().getGlobalNumEntries ();
 
  531   const Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>&
 
  541   typedef Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> MV;
 
  545   void allocateSolvers ();
 
  546   void allocate_L_and_U ();
 
  559   typedef Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> vec_type;
 
  592   mutable int numApply_;
 
  594   double initializeTime_;
 
  596   mutable double applyTime_;
 
  611   template<
class MatrixType, 
class NodeType>
 
  612   struct setLocalSolveParams{
 
  620 template <
class MatrixType>
 
  621 template <
typename NewMatrixType>
 
  630   typedef typename NewMatrixType::scalar_type new_scalar_type;
 
  631   typedef typename NewMatrixType::local_ordinal_type new_local_ordinal_type;
 
  632   typedef typename NewMatrixType::global_ordinal_type new_global_ordinal_type;
 
  633   typedef typename NewMatrixType::node_type new_node_type;
 
  634   typedef Tpetra::RowMatrix<new_scalar_type, new_local_ordinal_type,
 
  635     new_global_ordinal_type, new_node_type> new_row_matrix_type;
 
  638   RCP<new_riluk_type> new_riluk = 
rcp (
new new_riluk_type (A_newnode));
 
  640   RCP<ParameterList> plClone = Teuchos::parameterList ();
 
  641   plClone = detail::setLocalSolveParams<NewMatrixType, new_node_type>::setParams (plClone);
 
  647   new_riluk->LevelOfFill_ = LevelOfFill_;
 
  649   new_riluk->isAllocated_ = isAllocated_;
 
  650   new_riluk->isInitialized_ = isInitialized_;
 
  651   new_riluk->isComputed_ = isComputed_;
 
  653   new_riluk->numInitialize_ = numInitialize_;
 
  654   new_riluk->numCompute_ = numCompute_;
 
  655   new_riluk->numApply_ =  numApply_;
 
  657   new_riluk->RelaxValue_ = RelaxValue_;
 
  658   new_riluk->Athresh_ = Athresh_;
 
  659   new_riluk->Rthresh_ = Rthresh_;
 
Mix-in interface for preconditioners that can change their matrix after construction. 
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
 
Declaration and definition of IlukGraph. 
 
magnitude_type getRelaxValue() const 
Get RILU(k) relaxation parameter. 
Definition: Ifpack2_RILUK_decl.hpp:498
 
magnitude_type getAbsoluteThreshold() const 
Get absolute threshold value. 
Definition: Ifpack2_RILUK_decl.hpp:501
 
Ifpack2::ScalingType enumerable type. 
 
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const 
Returns the Tpetra::Map object associated with the range of this operator. 
Definition: Ifpack2_RILUK_def.hpp:222
 
Ifpack2::Preconditioner< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::scalar_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::local_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::global_ordinal_type, Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::node_type >::magnitude_type Teuchos::ScalarTraits< Tpetra::RowMatrix< MatrixType::scalar_type, MatrixType::local_ordinal_type, MatrixType::global_ordinal_type, MatrixType::node_type >::scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry. 
Definition: Ifpack2_Preconditioner.hpp:111
 
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType. 
Definition: Ifpack2_RILUK_decl.hpp:261
 
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k). 
Definition: Ifpack2_RILUK_decl.hpp:574
 
const Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > & getD() const 
Return the diagonal entries of the ILU factorization. 
Definition: Ifpack2_RILUK_def.hpp:158
 
void initialize()
Initialize by computing the symbolic incomplete factorization. 
Definition: Ifpack2_RILUK_def.hpp:387
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class. 
Definition: Ifpack2_RILUK_decl.hpp:273
 
double getApplyTime() const 
Total time in seconds taken by all successful apply() calls for this object. 
Definition: Ifpack2_RILUK_decl.hpp:371
 
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned. 
Definition: Ifpack2_RILUK_def.hpp:105
 
std::string description() const 
A one-line description of this object. 
Definition: Ifpack2_RILUK_def.hpp:974
 
void setParameters(const Teuchos::ParameterList ¶ms)
Definition: Ifpack2_RILUK_def.hpp:269
 
ILU(k) factorization of a given Tpetra::RowMatrix. 
Definition: Ifpack2_RILUK_decl.hpp:243
 
size_t getNodeSmootherComplexity() const 
Get a rough estimate of cost per iteration. 
Definition: Ifpack2_RILUK_def.hpp:185
 
Tpetra::CombineMode getOverlapMode()
Get overlap mode type. 
Definition: Ifpack2_RILUK_decl.hpp:510
 
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > getGraph() const 
Return the Ifpack2::IlukGraph associated with this factored matrix. 
Definition: Ifpack2_RILUK_decl.hpp:523
 
Teuchos::RCP< const crs_matrix_type > getCrsMatrix() const 
Return the input matrix A as a Tpetra::CrsMatrix, if possible; else throws. 
Definition: Ifpack2_RILUK_def.hpp:347
 
bool isComputed() const 
Whether compute() has been called on this object. 
Definition: Ifpack2_RILUK_decl.hpp:345
 
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class for representing L and U. 
Definition: Ifpack2_RILUK_decl.hpp:276
 
const crs_matrix_type & getU() const 
Return the U factor of the ILU factorization. 
Definition: Ifpack2_RILUK_def.hpp:172
 
MatrixType::node_type node_type
The Node type used by the input MatrixType. 
Definition: Ifpack2_RILUK_decl.hpp:264
 
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const 
Returns the Tpetra::Map object associated with the domain of this operator. 
Definition: Ifpack2_RILUK_def.hpp:203
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType. 
Definition: Ifpack2_RILUK_decl.hpp:258
 
virtual ~RILUK()
Destructor (declared virtual for memory safety). 
Definition: Ifpack2_RILUK_def.hpp:94
 
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > U_solver_
Sparse triangular solver for U. 
Definition: Ifpack2_RILUK_decl.hpp:580
 
int getLevelOfFill() const 
Get level of fill (the "k" in ILU(k)). 
Definition: Ifpack2_RILUK_decl.hpp:507
 
Teuchos::RCP< RILUK< NewMatrixType > > clone(const Teuchos::RCP< const NewMatrixType > &A_newnode) const 
Clone preconditioner to a new node type. 
Definition: Ifpack2_RILUK_decl.hpp:624
 
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > Graph_
The ILU(k) graph. 
Definition: Ifpack2_RILUK_decl.hpp:567
 
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:107
 
void compute()
Compute the (numeric) incomplete factorization. 
Definition: Ifpack2_RILUK_def.hpp:644
 
magnitude_type getRelativeThreshold() const 
Get relative threshold value. 
Definition: Ifpack2_RILUK_decl.hpp:504
 
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const 
Apply the (inverse of the) incomplete factorization to X, resulting in Y. 
Definition: Ifpack2_RILUK_def.hpp:826
 
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry. 
Definition: Ifpack2_RILUK_decl.hpp:267
 
Construct a level filled graph for use in computing an ILU(k) incomplete factorization. 
Definition: Ifpack2_IlukGraph.hpp:97
 
Declaration of interface for preconditioners that can change their matrix after construction. 
 
bool isInitialized() const 
Whether initialize() has been called on this object. 
Definition: Ifpack2_RILUK_decl.hpp:341
 
Teuchos::RCP< LocalSparseTriangularSolver< row_matrix_type > > L_solver_
Sparse triangular solver for L. 
Definition: Ifpack2_RILUK_decl.hpp:576
 
Teuchos::RCP< const row_matrix_type > A_
The (original) input matrix for which to compute ILU(k). 
Definition: Ifpack2_RILUK_decl.hpp:562
 
Teuchos::RCP< const row_matrix_type > getMatrix() const 
Get the input matrix. 
Definition: Ifpack2_RILUK_def.hpp:340
 
double getComputeTime() const 
Total time in seconds taken by all successful compute() calls for this object. 
Definition: Ifpack2_RILUK_decl.hpp:367
 
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization. 
Definition: Ifpack2_RILUK_decl.hpp:582
 
int getNumCompute() const 
Number of successful compute() calls for this object. 
Definition: Ifpack2_RILUK_decl.hpp:354
 
int getNumApply() const 
Number of successful apply() calls for this object. 
Definition: Ifpack2_RILUK_decl.hpp:358
 
Tpetra::global_size_t getGlobalNumEntries() const 
Returns the number of nonzero entries in the global graph. 
Definition: Ifpack2_RILUK_decl.hpp:518
 
double getInitializeTime() const 
Total time in seconds taken by all successful initialize() calls for this object. ...
Definition: Ifpack2_RILUK_decl.hpp:363
 
const crs_matrix_type & getL() const 
Return the L factor of the ILU factorization. 
Definition: Ifpack2_RILUK_def.hpp:141
 
Teuchos::RCP< const row_matrix_type > A_local_
The matrix whos numbers are used to to compute ILU(k). The graph may be computed using a crs_matrix_t...
Definition: Ifpack2_RILUK_decl.hpp:571
 
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k). 
Definition: Ifpack2_RILUK_decl.hpp:578
 
int getNumInitialize() const 
Number of successful initialize() calls for this object. 
Definition: Ifpack2_RILUK_decl.hpp:350
 
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType. 
Definition: Ifpack2_RILUK_decl.hpp:255