43#ifndef _IFPACK_CRSRILUK_H_
44#define _IFPACK_CRSRILUK_H_
46#if defined(Ifpack_SHOW_DEPRECATED_WARNINGS)
48#warning "The Ifpack package is deprecated"
52#include "Ifpack_ConfigDefs.h"
54#include "Ifpack_IlukGraph.h"
55#include "Epetra_ConfigDefs.h"
56#include "Epetra_CompObject.h"
57#include "Epetra_Operator.h"
58#include "Epetra_CrsMatrix.h"
59#include "Epetra_Object.h"
60#include "Epetra_MultiVector.h"
61#include "Epetra_Vector.h"
62#include "Epetra_Map.h"
67#include "Teuchos_RefCountPtr.hpp"
253#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
280 int SetParameters(
const Teuchos::ParameterList& parameterlist,
281 bool cerr_warning_if_unused=
false);
336 int Condest(
bool Trans,
double & ConditionNumberEstimate)
const;
352#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
366 long long NumGlobalRows64()
const {
return(
Graph().NumGlobalRows64());};
367 long long NumGlobalCols64()
const {
return(
Graph().NumGlobalCols64());};
368 long long NumGlobalNonzeros64()
const {
return(
L().NumGlobalNonzeros64()+
U().NumGlobalNonzeros64());};
369 virtual long long NumGlobalBlockDiagonals64()
const {
return(
Graph().NumGlobalBlockDiagonals64());};
387#ifndef EPETRA_NO_32BIT_GLOBAL_INDICES
390 long long IndexBase64()
const {
return(
Graph().IndexBase64());};
419 int SetUseTranspose(
bool UseTranspose_in) {UseTranspose_ = UseTranspose_in;
return(0);};
473 void SetFactored(
bool Flag) {Factored_ = Flag;};
474 void SetValuesInitialized(
bool Flag) {ValuesInitialized_ = Flag;};
475 bool Allocated()
const {
return(Allocated_);};
476 int SetAllocated(
bool Flag) {Allocated_ = Flag;
return(0);};
485 int BlockMap2PointMap(
const Epetra_BlockMap & BlockMap, Teuchos::RefCountPtr<Epetra_Map>* PointMap);
486 int GenerateXY(
bool Trans,
488 Teuchos::RefCountPtr<Epetra_MultiVector>* Xout,
489 Teuchos::RefCountPtr<Epetra_MultiVector>* Yout)
const;
490 bool UserMatrixIsVbr_;
491 bool UserMatrixIsCrs_;
494 Teuchos::RefCountPtr<Epetra_Map> IlukRowMap_;
495 Teuchos::RefCountPtr<Epetra_Map> IlukDomainMap_;
496 Teuchos::RefCountPtr<Epetra_Map> IlukRangeMap_;
497 Teuchos::RefCountPtr<const Epetra_Map> U_DomainMap_;
498 Teuchos::RefCountPtr<const Epetra_Map> L_RangeMap_;
500 Teuchos::RefCountPtr<Epetra_CrsMatrix> L_;
501 Teuchos::RefCountPtr<Epetra_CrsMatrix> U_;
502 Teuchos::RefCountPtr<Epetra_CrsGraph> L_Graph_;
503 Teuchos::RefCountPtr<Epetra_CrsGraph> U_Graph_;
504 Teuchos::RefCountPtr<Epetra_Vector> D_;
509 bool ValuesInitialized_;
514 mutable double Condest_;
516 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlapX_;
517 mutable Teuchos::RefCountPtr<Epetra_MultiVector> OverlapY_;
518 mutable Teuchos::RefCountPtr<Epetra_MultiVector> VbrX_;
519 mutable Teuchos::RefCountPtr<Epetra_MultiVector> VbrY_;
Ifpack_ScalingType enumerable type.
virtual const char * Label() const
Ifpack_CrsRiluk: A class for constructing and using an incomplete lower/upper (ILU) factorization of ...
const Epetra_Map & OperatorDomainMap() const
Returns the Epetra_Map object associated with the domain of this operator.
const Epetra_Vector & D() const
Returns the address of the D factor associated with this factored matrix.
virtual int NumMyBlockDiagonals() const
Returns the number of diagonal entries found in the local input graph.
int NumGlobalCols() const
Returns the number of global matrix columns.
virtual ~Ifpack_CrsRiluk()
Ifpack_CrsRiluk Destructor.
int Multiply(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of multiplying U, D and L in that order on an Epetra_MultiVector X in Y.
int Factor()
Compute ILU factors L and U using the specified graph, diagonal perturbation thresholds and relaxatio...
int SetParameters(const Teuchos::ParameterList ¶meterlist, bool cerr_warning_if_unused=false)
Set parameters using a Teuchos::ParameterList object.
void SetRelaxValue(double RelaxValue)
Set RILU(k) relaxation parameter.
int NumMyRows() const
Returns the number of local matrix rows.
friend std::ostream & operator<<(std::ostream &os, const Ifpack_CrsRiluk &A)
<< operator will work for Ifpack_CrsRiluk.
int SetUseTranspose(bool UseTranspose_in)
If set true, transpose of this operator will be applied.
const Epetra_CrsMatrix & U() const
Returns the address of the L factor associated with this factored matrix.
Epetra_CombineMode GetOverlapMode()
Get overlap mode type.
bool Factored() const
If factor is completed, this query returns true, otherwise it returns false.
int IndexBase() const
Returns the index base for row and column indices for this graph.
bool UseTranspose() const
Returns the current UseTranspose setting.
bool HasNormInf() const
Returns false because this class cannot compute an Inf-norm.
void SetAbsoluteThreshold(double Athresh)
Set absolute threshold value.
int InitValues(const Epetra_CrsMatrix &A)
Initialize L and U with values from user matrix A.
double GetRelaxValue()
Get RILU(k) relaxation parameter.
const Epetra_Map & OperatorRangeMap() const
Returns the Epetra_Map object associated with the range of this operator.
const char * Label() const
Returns a character string describing the operator.
int NumGlobalRows() const
Returns the number of global matrix rows.
int NumGlobalNonzeros() const
Returns the number of nonzero entries in the global graph.
virtual int NumGlobalBlockDiagonals() const
Returns the number of diagonal entries found in the global input graph.
virtual int NumMyDiagonals() const
Returns the number of nonzero diagonal values found in matrix.
int ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
void SetRelativeThreshold(double Rthresh)
Set relative threshold value.
const Epetra_Comm & Comm() const
Returns the Epetra_BlockMap object associated with the range of this matrix operator.
int Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
bool ValuesInitialized() const
If values have been initialized, this query returns true, otherwise it returns false.
double GetAbsoluteThreshold()
Get absolute threshold value.
const Epetra_CrsMatrix & L() const
Returns the address of the L factor associated with this factored matrix.
int NumMyNonzeros() const
Returns the number of nonzero entries in the local graph.
double GetRelativeThreshold()
Get relative threshold value.
int NumMyCols() const
Returns the number of local matrix columns.
void SetOverlapMode(Epetra_CombineMode OverlapMode)
Set overlap mode type.
int Condest(bool Trans, double &ConditionNumberEstimate) const
Returns the maximum over all the condition number estimate for each local ILU set of factors.
int Solve(bool Trans, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Returns the result of a Ifpack_CrsRiluk forward/back solve on a Epetra_MultiVector X in Y (works for ...
double NormInf() const
Returns 0.0 because this class cannot compute Inf-norm.
const Ifpack_IlukGraph & Graph() const
returns the address of the Ifpack_IlukGraph associated with this factored matrix.
Ifpack_IlukGraph: A class for constructing level filled graphs for use with ILU(k) class precondition...