PyTrilinos  Development
 All Classes Functions
Public Member Functions | Public Attributes
PyTrilinos.IFPACK.ILU Class Reference
Inheritance diagram for PyTrilinos.IFPACK.ILU:
Inheritance graph
[legend]
Collaboration diagram for PyTrilinos.IFPACK.ILU:
Collaboration graph
[legend]

List of all members.

Public Member Functions

def __init__
def Initialize
def IsInitialized
def Compute
def IsComputed
def SetParameters
def SetUseTranspose
def Apply
def Multiply
def ApplyInverse
def Condest
def L
def D
def U
def Label
def SetLabel
def NormInf
def HasNormInf
def UseTranspose
def OperatorDomainMap
def OperatorRangeMap
def Comm
def Matrix
def NumInitialize
def NumCompute
def NumApplyInverse
def InitializeTime
def ComputeTime
def ApplyInverseTime
def InitializeFlops
def ComputeFlops
def ApplyInverseFlops

Public Attributes

 this

Detailed Description

Ifpack_ILU: A class for constructing and using an incomplete
lower/upper (ILU) factorization of a given Epetra_RowMatrix.

The Ifpack_ILU class computes a "Relaxed" ILU factorization with
level k fill of a given Epetra_RowMatrix.

Please refer to ifp_ilu for a general description of the ILU
algorithm.

The complete list of supported parameters is reported in page
ifp_params.

Mike Heroux, Marzio Sala, SNL 9214.

C++ includes: Ifpack_ILU.h 

Constructor & Destructor Documentation

def PyTrilinos.IFPACK.ILU.__init__ (   self,
  args 
)
__init__(Ifpack_ILU self, RowMatrix A) -> ILU

Ifpack_ILU::Ifpack_ILU(Epetra_RowMatrix *A)

Constructor. 

Reimplemented from PyTrilinos.Epetra.Operator.


Member Function Documentation

def PyTrilinos.IFPACK.ILU.Apply (   self,
  args 
)
Apply(ILU self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int

int Ifpack_ILU::Apply(const
Epetra_MultiVector &X, Epetra_MultiVector &Y) const 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.ApplyInverse (   self,
  args 
)
ApplyInverse(ILU self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int

int
Ifpack_ILU::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.

In this implementation, we use several existing attributes to
determine how virtual method ApplyInverse() should call the concrete
method Solve(). We pass in the UpperTriangular(), the
Epetra_CrsMatrix::UseTranspose(), and NoDiagonal() methods. The most
notable warning is that if a matrix has no diagonal values we assume
that there is an implicit unit diagonal that should be accounted for
when doing a triangular solve.

Parameters:
-----------

X:  - (In) A Epetra_MultiVector of dimension NumVectors to solve for.

Out:  Y - (Out) A Epetra_MultiVector of dimension NumVectors
containing result.

Integer error code, set to 0 if successful. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.ApplyInverseFlops (   self,
  args 
)
ApplyInverseFlops(ILU self) -> double

virtual double
Ifpack_ILU::ApplyInverseFlops() const

Returns the number of flops in the application of the preconditioner.

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.ApplyInverseTime (   self,
  args 
)
ApplyInverseTime(ILU self) -> double

virtual double
Ifpack_ILU::ApplyInverseTime() const

Returns the time spent in ApplyInverse(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.Comm (   self,
  args 
)
Comm(ILU self) -> Comm

const Epetra_Comm&
Ifpack_ILU::Comm() const

Returns the Epetra_BlockMap object associated with the range of this
matrix operator. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.Compute (   self,
  args 
)
Compute(ILU self) -> int

int Ifpack_ILU::Compute()

Compute ILU factors L and U using the specified graph, diagonal
perturbation thresholds and relaxation parameters.

This function computes the ILU(k) factors L and U using the current:
Ifpack_IlukGraph specifying the structure of L and U.

Value for the ILU(k) relaxation parameter.

Value for the a priori diagonal threshold values.  InitValues() must
be called before the factorization can proceed. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.ComputeFlops (   self,
  args 
)
ComputeFlops(ILU self) -> double

virtual double
Ifpack_ILU::ComputeFlops() const

Returns the number of flops in the computation phase. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.ComputeTime (   self,
  args 
)
ComputeTime(ILU self) -> double

virtual double
Ifpack_ILU::ComputeTime() const

Returns the time spent in Compute(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.Condest (   self,
  args 
)
Condest(ILU self, Ifpack_CondestType const CT=Ifpack_Cheap, int const MaxIters=1550, double const Tol=1e-9, 
    RowMatrix Matrix_in=None) -> double

double
Ifpack_ILU::Condest() const

Returns the computed estimated condition number, or -1.0 if not
computed. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.D (   self,
  args 
)
D(ILU self) -> Epetra_Vector

const Epetra_Vector&
Ifpack_ILU::D() const

Returns the address of the D factor associated with this factored
matrix. 
def PyTrilinos.IFPACK.ILU.HasNormInf (   self,
  args 
)
HasNormInf(ILU self) -> bool

bool
Ifpack_ILU::HasNormInf() const

Returns false because this class cannot compute an Inf-norm. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.Initialize (   self,
  args 
)
Initialize(ILU self) -> int

int
Ifpack_ILU::Initialize()

Initialize the preconditioner, does not touch matrix values. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.InitializeFlops (   self,
  args 
)
InitializeFlops(ILU self) -> double

virtual double
Ifpack_ILU::InitializeFlops() const

Returns the number of flops in the initialization phase. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.InitializeTime (   self,
  args 
)
InitializeTime(ILU self) -> double

virtual double
Ifpack_ILU::InitializeTime() const

Returns the time spent in Initialize(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.IsComputed (   self,
  args 
)
IsComputed(ILU self) -> bool

bool
Ifpack_ILU::IsComputed() const

If factor is completed, this query returns true, otherwise it returns
false. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.IsInitialized (   self,
  args 
)
IsInitialized(ILU self) -> bool

bool
Ifpack_ILU::IsInitialized() const

Returns true if the preconditioner has been successfully initialized.

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.L (   self,
  args 
)
L(ILU self) -> CrsMatrix

const Epetra_CrsMatrix&
Ifpack_ILU::L() const

Returns the address of the L factor associated with this factored
matrix. 
def PyTrilinos.IFPACK.ILU.Label (   self,
  args 
)
Label(ILU self) -> char const *

const char*
Ifpack_ILU::Label() const

Returns a character string describing the operator. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.Matrix (   self,
  args 
)
Matrix(ILU self) -> RowMatrix

const Epetra_RowMatrix&
Ifpack_ILU::Matrix() const

Returns a reference to the matrix to be preconditioned. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.Multiply (   self,
  args 
)
Multiply(ILU self, bool Trans, Epetra_MultiVector X, Epetra_MultiVector Y) -> int

int
Ifpack_ILU::Multiply(bool Trans, const Epetra_MultiVector &X,
Epetra_MultiVector &Y) const 
def PyTrilinos.IFPACK.ILU.NormInf (   self,
  args 
)
NormInf(ILU self) -> double

double
Ifpack_ILU::NormInf() const

Returns 0.0 because this class cannot compute Inf-norm. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.NumApplyInverse (   self,
  args 
)
NumApplyInverse(ILU self) -> int

virtual int
Ifpack_ILU::NumApplyInverse() const

Returns the number of calls to ApplyInverse(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.NumCompute (   self,
  args 
)
NumCompute(ILU self) -> int

virtual int
Ifpack_ILU::NumCompute() const

Returns the number of calls to Compute(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.NumInitialize (   self,
  args 
)
NumInitialize(ILU self) -> int

virtual int
Ifpack_ILU::NumInitialize() const

Returns the number of calls to Initialize(). 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.OperatorDomainMap (   self,
  args 
)
OperatorDomainMap(ILU self) -> Map

const
Epetra_Map& Ifpack_ILU::OperatorDomainMap() const

Returns the Epetra_Map object associated with the domain of this
operator. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.OperatorRangeMap (   self,
  args 
)
OperatorRangeMap(ILU self) -> Map

const
Epetra_Map& Ifpack_ILU::OperatorRangeMap() const

Returns the Epetra_Map object associated with the range of this
operator. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.SetLabel (   self,
  args 
)
SetLabel(ILU self, char const * Label_in) -> int

int
Ifpack_ILU::SetLabel(const char *Label_in)

Sets label for this object. 
def PyTrilinos.IFPACK.ILU.SetParameters (   self,
  args 
)
SetParameters(ILU self, ParameterList parameterlist) -> int

int
Ifpack_ILU::SetParameters(Teuchos::ParameterList &parameterlist)

Set parameters using a Teuchos::ParameterList object. 

Reimplemented from PyTrilinos.IFPACK.Preconditioner.

def PyTrilinos.IFPACK.ILU.SetUseTranspose (   self,
  args 
)
SetUseTranspose(ILU self, bool UseTranspose_in) -> int

int
Ifpack_ILU::SetUseTranspose(bool UseTranspose_in)

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used
implicitly. Setting this flag affects only the Apply() and
ApplyInverse() methods. If the implementation of this interface does
not support transpose use, this method should return a value of -1.

Parameters:
-----------

UseTranspose_in:  - (In) If true, multiply by the transpose of
operator, otherwise just use operator.

Always returns 0. 

Reimplemented from PyTrilinos.Epetra.Operator.

def PyTrilinos.IFPACK.ILU.U (   self,
  args 
)
U(ILU self) -> CrsMatrix

const Epetra_CrsMatrix&
Ifpack_ILU::U() const

Returns the address of the L factor associated with this factored
matrix. 
def PyTrilinos.IFPACK.ILU.UseTranspose (   self,
  args 
)
UseTranspose(ILU self) -> bool

bool
Ifpack_ILU::UseTranspose() const

Returns the current UseTranspose setting. 

Reimplemented from PyTrilinos.Epetra.Operator.


The documentation for this class was generated from the following file:
 All Classes Functions