

Public Member Functions | |
| def | __init__ |
| def | SetAbsoluteThreshold |
| def | SetRelativeThreshold |
| def | SetParameters |
| def | SetParameter |
| def | Matrix |
| def | IsInitialized |
| def | Initialize |
| def | Compute |
| def | ComputeSetup |
| def | IsComputed |
| def | ApplyInverse |
| def | Apply |
| def | Condest |
| def | GetAbsoluteThreshold |
| def | GetRelativeThreshold |
| def | NumGlobalNonzeros |
| def | NumMyNonzeros |
| def | D |
| def | U |
| def | SetUseTranspose |
| def | NormInf |
| def | HasNormInf |
| def | UseTranspose |
| def | OperatorDomainMap |
| def | OperatorRangeMap |
| def | Comm |
| def | Label |
| def | SetLabel |
| def | NumInitialize |
| def | NumCompute |
| def | NumApplyInverse |
| def | InitializeTime |
| def | ComputeTime |
| def | ApplyInverseTime |
| def | InitializeFlops |
| def | ComputeFlops |
| def | ApplyInverseFlops |
Public Attributes | |
| this | |
Ifpack_IC: A class for constructing and using an incomplete Cholesky factorization of a given Epetra_RowMatrix. The Ifpack_IC class computes a threshold based incomplete LDL^T factorization of a given Epetra_RowMatrix. The factorization that is produced is a function of several parameters: Maximum number of entries per row/column in factor - The factorization will contain at most this number of nonzero terms in each row/column of the factorization. Diagonal perturbation - Prior to computing the factorization, it is possible to modify the diagonal entries of the matrix for which the factorization will be computing. If the absolute and relative perturbation values are zero and one, respectively, the factorization will be compute for the original user matrix A. Otherwise, the factorization will computed for a matrix that differs from the original user matrix in the diagonal values only. Details can be found in ifp_diag_pert. C++ includes: Ifpack_IC.h
| def PyTrilinos.IFPACK.IC.__init__ | ( | self, | |
| args | |||
| ) |
__init__(Ifpack_IC self, RowMatrix A) -> IC Ifpack_IC::Ifpack_IC(Epetra_RowMatrix *A) Ifpack_IC constuctor with variable number of indices per row. Creates a Ifpack_IC object and allocates storage. Parameters: ----------- In: A - User matrix to be factored. In: Graph - Graph generated by Ifpack_IlukGraph.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.Apply | ( | self, | |
| args | |||
| ) |
Apply(IC self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int Ifpack_IC::Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.ApplyInverse | ( | self, | |
| args | |||
| ) |
ApplyInverse(IC self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int Ifpack_IC::ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const Returns the result of a Ifpack_IC forward/back solve on a Epetra_MultiVector X in Y. Parameters: ----------- In: Trans -If true, solve transpose problem. In: X - A Epetra_MultiVector of dimension NumVectors to solve for. Out: Y -A Epetra_MultiVector of dimension NumVectorscontaining result. Integer error code, set to 0 if successful.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.ApplyInverseFlops | ( | self, | |
| args | |||
| ) |
ApplyInverseFlops(IC self) -> double virtual double Ifpack_IC::ApplyInverseFlops() const Returns the number of flops in the application of the preconditioner.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.ApplyInverseTime | ( | self, | |
| args | |||
| ) |
ApplyInverseTime(IC self) -> double virtual double Ifpack_IC::ApplyInverseTime() const Returns the time spent in ApplyInverse().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.Comm | ( | self, | |
| args | |||
| ) |
Comm(IC self) -> Comm const Epetra_Comm& Ifpack_IC::Comm() const Returns the Epetra_BlockMap object associated with the range of this matrix operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.Compute | ( | self, | |
| args | |||
| ) |
Compute(IC self) -> int int Ifpack_IC::Compute() Compute IC factor U using the specified graph, diagonal perturbation thresholds and relaxation parameters. This function computes the RILU(k) factors L and U using the current: Ifpack_IlukGraph specifying the structure of L and U. Value for the RILU(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.IC.ComputeFlops | ( | self, | |
| args | |||
| ) |
ComputeFlops(IC self) -> double virtual double Ifpack_IC::ComputeFlops() const Returns the number of flops in the computation phase.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.ComputeSetup | ( | self, | |
| args | |||
| ) |
ComputeSetup(IC self) -> int int Ifpack_IC::ComputeSetup()
| def PyTrilinos.IFPACK.IC.ComputeTime | ( | self, | |
| args | |||
| ) |
ComputeTime(IC self) -> double virtual double Ifpack_IC::ComputeTime() const Returns the time spent in Compute().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.Condest | ( | self, | |
| args | |||
| ) |
Condest(IC self, Ifpack_CondestType const CT=Ifpack_Cheap, int const MaxIters=1550, double const Tol=1e-9,
RowMatrix Matrix_in=None) -> double
double
Ifpack_IC::Condest() const
Returns the computed condition number estimate, or -1.0 if not
computed.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.D | ( | self, | |
| args | |||
| ) |
D(IC self) -> Epetra_Vector const Epetra_Vector& Ifpack_IC::D() const Returns the address of the D factor associated with this factored matrix.
| def PyTrilinos.IFPACK.IC.GetAbsoluteThreshold | ( | self, | |
| args | |||
| ) |
GetAbsoluteThreshold(IC self) -> double double Ifpack_IC::GetAbsoluteThreshold() Get absolute threshold value.
| def PyTrilinos.IFPACK.IC.GetRelativeThreshold | ( | self, | |
| args | |||
| ) |
GetRelativeThreshold(IC self) -> double double Ifpack_IC::GetRelativeThreshold() Get relative threshold value.
| def PyTrilinos.IFPACK.IC.HasNormInf | ( | self, | |
| args | |||
| ) |
HasNormInf(IC self) -> bool bool Ifpack_IC::HasNormInf() const Returns false because this class cannot compute an Inf-norm.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.Initialize | ( | self, | |
| args | |||
| ) |
Initialize(IC self) -> int int Ifpack_IC::Initialize() Initialize L and U with values from user matrix A. Copies values from the user's matrix into the nonzero pattern of L and U. Parameters: ----------- In: A - User matrix to be factored. WARNING: The graph of A must be identical to the graph passed in to Ifpack_IlukGraph constructor.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.InitializeFlops | ( | self, | |
| args | |||
| ) |
InitializeFlops(IC self) -> double virtual double Ifpack_IC::InitializeFlops() const Returns the number of flops in the initialization phase.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.InitializeTime | ( | self, | |
| args | |||
| ) |
InitializeTime(IC self) -> double virtual double Ifpack_IC::InitializeTime() const Returns the time spent in Initialize().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.IsComputed | ( | self, | |
| args | |||
| ) |
IsComputed(IC self) -> bool bool Ifpack_IC::IsComputed() const If factor is completed, this query returns true, otherwise it returns false.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.IsInitialized | ( | self, | |
| args | |||
| ) |
IsInitialized(IC self) -> bool bool Ifpack_IC::IsInitialized() const Returns true if the preconditioner has been successfully initialized, false otherwise.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.Label | ( | self, | |
| args | |||
| ) |
Label(IC self) -> char const * const char* Ifpack_IC::Label() const
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.Matrix | ( | self, | |
| args | |||
| ) |
Matrix(IC self) -> RowMatrix Matrix(IC self) -> RowMatrix Epetra_RowMatrix& Ifpack_IC::Matrix()
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.NormInf | ( | self, | |
| args | |||
| ) |
NormInf(IC self) -> double double Ifpack_IC::NormInf() const Returns 0.0 because this class cannot compute Inf-norm.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.NumApplyInverse | ( | self, | |
| args | |||
| ) |
NumApplyInverse(IC self) -> int virtual int Ifpack_IC::NumApplyInverse() const Returns the number of calls to ApplyInverse().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.NumCompute | ( | self, | |
| args | |||
| ) |
NumCompute(IC self) -> int virtual int Ifpack_IC::NumCompute() const Returns the number of calls to Compute().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.NumGlobalNonzeros | ( | self, | |
| args | |||
| ) |
NumGlobalNonzeros(IC self) -> int int Ifpack_IC::NumGlobalNonzeros() const Returns the number of nonzero entries in the global graph.
| def PyTrilinos.IFPACK.IC.NumInitialize | ( | self, | |
| args | |||
| ) |
NumInitialize(IC self) -> int virtual int Ifpack_IC::NumInitialize() const Returns the number of calls to Initialize().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.NumMyNonzeros | ( | self, | |
| args | |||
| ) |
NumMyNonzeros(IC self) -> int int Ifpack_IC::NumMyNonzeros() const Returns the number of nonzero entries in the local graph.
| def PyTrilinos.IFPACK.IC.OperatorDomainMap | ( | self, | |
| args | |||
| ) |
OperatorDomainMap(IC self) -> Map const Epetra_Map& Ifpack_IC::OperatorDomainMap() const Returns the Epetra_Map object associated with the domain of this operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.OperatorRangeMap | ( | self, | |
| args | |||
| ) |
OperatorRangeMap(IC self) -> Map const Epetra_Map& Ifpack_IC::OperatorRangeMap() const Returns the Epetra_Map object associated with the range of this operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.SetAbsoluteThreshold | ( | self, | |
| args | |||
| ) |
SetAbsoluteThreshold(IC self, double Athresh) void Ifpack_IC::SetAbsoluteThreshold(double Athresh) Set absolute threshold value.
| def PyTrilinos.IFPACK.IC.SetLabel | ( | self, | |
| args | |||
| ) |
SetLabel(IC self, char const * Label_in) -> int int Ifpack_IC::SetLabel(const char *Label_in)
| def PyTrilinos.IFPACK.IC.SetParameter | ( | self, | |
| args | |||
| ) |
SetParameter(IC self, string const Name, int const Value) -> int SetParameter(IC self, string const Name, double const Value) -> int int Ifpack_IC::SetParameter(const string Name, const double Value)
| def PyTrilinos.IFPACK.IC.SetParameters | ( | self, | |
| args | |||
| ) |
SetParameters(IC self, ParameterList parameterlis) -> int int Ifpack_IC::SetParameters(Teuchos::ParameterList ¶meterlis) Set parameters using a Teuchos::ParameterList object.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.IC.SetRelativeThreshold | ( | self, | |
| args | |||
| ) |
SetRelativeThreshold(IC self, double Rthresh) void Ifpack_IC::SetRelativeThreshold(double Rthresh) Set relative threshold value.
| def PyTrilinos.IFPACK.IC.SetUseTranspose | ( | self, | |
| args | |||
| ) |
SetUseTranspose(IC self, bool UseTranspose_in) -> int int Ifpack_IC::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: ----------- In: UseTranspose_in -If true, multiply by the transpose of operator, otherwise just use operator. Always returns 0.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.IC.U | ( | self, | |
| args | |||
| ) |
U(IC self) -> CrsMatrix const Epetra_CrsMatrix& Ifpack_IC::U() const Returns the address of the U factor associated with this factored matrix.
| def PyTrilinos.IFPACK.IC.UseTranspose | ( | self, | |
| args | |||
| ) |
UseTranspose(IC self) -> bool bool Ifpack_IC::UseTranspose() const Returns the current UseTranspose setting.
Reimplemented from PyTrilinos.Epetra.Operator.
1.7.6.1