

Public Member Functions | |
| def | __init__ |
| def | SetUseTranspose |
| def | Apply |
| def | ApplyInverse |
| def | NormInf |
| def | Label |
| def | UseTranspose |
| def | HasNormInf |
| def | Comm |
| def | OperatorDomainMap |
| def | OperatorRangeMap |
| def | Initialize |
| def | IsInitialized |
| def | IsComputed |
| def | Compute |
| def | Matrix |
| def | Condest |
| def | SetParameters |
| def | NumInitialize |
| def | NumCompute |
| def | NumApplyInverse |
| def | InitializeTime |
| def | ComputeTime |
| def | ApplyInverseTime |
| def | InitializeFlops |
| def | ComputeFlops |
| def | ApplyInverseFlops |
Public Attributes | |
| this | |
Ifpack_PointRelaxation: a class to define point relaxation
preconditioners of for Epetra_RowMatrix's.
The Ifpack_PointRelaxation class enables the construction of point
relaxation preconditioners of an Epetra_RowMatrix.
Ifpack_PointRelaxation is derived from the Ifpack_Preconditioner
class, which is itself derived from Epetra_Operator. Therefore this
object can be used as preconditioner everywhere an ApplyInverse()
method is required in the preconditioning step.
This class enables the construction of the following simple
preconditioners: Jacobi;
Gauss-Seidel;
symmetric Gauss-Seidel.
We now briefly describe the main features of the above
preconditioners. Consider a linear system of type \\[ A x = b, \\]
where $A$ is a square, real matrix, and $x, b$ are two real vectors.
We begin with the decomposition \\[ A = D - E - F \\] where $D$ is
the diagonal of A, $-E$ is the strict lower part, and $-F$ is the
strict upper part. It is assumed that the diagonal entries of $A$ are
different from zero.
Given an starting solution $x_0$, an iteration of the (damped) Jacobi
method can be written in matrix form as follows: \\[ x_{k+1} =
\\omega D^{-1}(E + F) x_k + D_{-1}b, \\] for $k < k_{max}$, and
$\\omega $ a damping parameter.
Using Ifpack_Jacobi, the user can apply the specified number of sweeps
( $k_{max}$), and the damping parameter. If only one sweep is used,
then the class simply applies the inverse of the diagonal of A to the
input vector.
Given an starting solution $x_0$, an iteration of the (damped)
GaussSeidel method can be written in matrix form as follows: \\[ (D
- E) x_{k+1} = \\omega F x_k + b, \\] for $k < k_{max}$, and
$\\omega $ a damping parameter. Equivalently, the Gauss-Seidel
preconditioner can be defined as \\[ P_{GS}^{-1} = (D - E)^{-1}.
\\] Clearly, the role of E and F can be interchanged. However,
Ifpack_GaussSeidel does not consider backward Gauss-Seidel methods.
For a list of supported parameters, please refer to page ifp_params.
The complete list of supported parameters is reported in page
ifp_params. For a presentation of basic relaxation schemes, please
refer to page Ifpack_PointRelaxation.
Marzio Sala, SNL 9214.
C++ includes: Ifpack_PointRelaxation.h
| def PyTrilinos.IFPACK.PointRelaxation.__init__ | ( | self, | |
| args | |||
| ) |
__init__(Ifpack_PointRelaxation self, RowMatrix Matrix) -> PointRelaxation Ifpack_PointRelaxation::Ifpack_PointRelaxation(const Epetra_RowMatrix *Matrix) Ifpack_PointRelaxation constructor with given Epetra_RowMatrix. Creates an instance of Ifpack_PointRelaxation class. Parameters: ----------- Matrix: - (In) Pointer to matrix to precondition.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.Apply | ( | self, | |
| args | |||
| ) |
Apply(PointRelaxation self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int virtual int Ifpack_PointRelaxation::Apply(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const Applies the matrix to an Epetra_MultiVector. Parameters: ----------- X: - (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix. Y: - (Out) A Epetra_MultiVector of dimension NumVectors containing the result. Integer error code, set to 0 if successful.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.ApplyInverse | ( | self, | |
| args | |||
| ) |
ApplyInverse(PointRelaxation self, Epetra_MultiVector X, Epetra_MultiVector Y) -> int int Ifpack_PointRelaxation::ApplyInverse(const Epetra_MultiVector &X, Epetra_MultiVector &Y) const Applies the preconditioner to X, returns the result in Y. Parameters: ----------- X: - (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned. Y: - (InOut) A Epetra_MultiVector of dimension NumVectors containing result. Integer error code, set to 0 if successful. WARNING: This routine is NOT AztecOO complaint.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.ApplyInverseFlops | ( | self, | |
| args | |||
| ) |
ApplyInverseFlops(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::ApplyInverseFlops() const Returns the number of flops for the application of the preconditioner.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.ApplyInverseTime | ( | self, | |
| args | |||
| ) |
ApplyInverseTime(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::ApplyInverseTime() const Returns the time spent in ApplyInverse().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.Comm | ( | self, | |
| args | |||
| ) |
Comm(PointRelaxation self) -> Comm const Epetra_Comm & Ifpack_PointRelaxation::Comm() const Returns a pointer to the Epetra_Comm communicator associated with this operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.Compute | ( | self, | |
| args | |||
| ) |
Compute(PointRelaxation self) -> int int Ifpack_PointRelaxation::Compute() Computes the preconditioners.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.ComputeFlops | ( | self, | |
| args | |||
| ) |
ComputeFlops(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::ComputeFlops() const Returns the number of flops in the computation phase.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.ComputeTime | ( | self, | |
| args | |||
| ) |
ComputeTime(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::ComputeTime() const Returns the time spent in Compute().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.Condest | ( | self, | |
| args | |||
| ) |
Condest(PointRelaxation self, Ifpack_CondestType const CT=Ifpack_Cheap, int const MaxIters=1550, double const Tol=1e-9,
RowMatrix Matrix=None) -> double
virtual
double Ifpack_PointRelaxation::Condest() const
Returns the condition number estimate, or -1.0 if not computed.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.HasNormInf | ( | self, | |
| args | |||
| ) |
HasNormInf(PointRelaxation self) -> bool virtual bool Ifpack_PointRelaxation::HasNormInf() const Returns true if the this object can provide an approximate Inf-norm, false otherwise.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.Initialize | ( | self, | |
| args | |||
| ) |
Initialize(PointRelaxation self) -> int int Ifpack_PointRelaxation::Initialize() Computes all it is necessary to initialize the preconditioner.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.InitializeFlops | ( | self, | |
| args | |||
| ) |
InitializeFlops(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::InitializeFlops() const Returns the number of flops in the initialization phase.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.InitializeTime | ( | self, | |
| args | |||
| ) |
InitializeTime(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::InitializeTime() const Returns the time spent in Initialize().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.IsComputed | ( | self, | |
| args | |||
| ) |
IsComputed(PointRelaxation self) -> bool virtual bool Ifpack_PointRelaxation::IsComputed() const Returns true if the preconditioner has been successfully computed.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.IsInitialized | ( | self, | |
| args | |||
| ) |
IsInitialized(PointRelaxation self) -> bool virtual bool Ifpack_PointRelaxation::IsInitialized() const Returns true if the preconditioner has been successfully initialized, false otherwise.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.Label | ( | self, | |
| args | |||
| ) |
Label(PointRelaxation self) -> char const * virtual const char* Ifpack_PointRelaxation::Label() const
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.Matrix | ( | self, | |
| args | |||
| ) |
Matrix(PointRelaxation self) -> RowMatrix virtual const Epetra_RowMatrix& Ifpack_PointRelaxation::Matrix() const Returns a pointer to the matrix to be preconditioned.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.NormInf | ( | self, | |
| args | |||
| ) |
NormInf(PointRelaxation self) -> double virtual double Ifpack_PointRelaxation::NormInf() const Returns the infinity norm of the global matrix (not implemented)
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.NumApplyInverse | ( | self, | |
| args | |||
| ) |
NumApplyInverse(PointRelaxation self) -> int virtual int Ifpack_PointRelaxation::NumApplyInverse() const Returns the number of calls to ApplyInverse().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.NumCompute | ( | self, | |
| args | |||
| ) |
NumCompute(PointRelaxation self) -> int virtual int Ifpack_PointRelaxation::NumCompute() const Returns the number of calls to Compute().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.NumInitialize | ( | self, | |
| args | |||
| ) |
NumInitialize(PointRelaxation self) -> int virtual int Ifpack_PointRelaxation::NumInitialize() const Returns the number of calls to Initialize().
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.OperatorDomainMap | ( | self, | |
| args | |||
| ) |
OperatorDomainMap(PointRelaxation self) -> Map const Epetra_Map & Ifpack_PointRelaxation::OperatorDomainMap() const Returns the Epetra_Map object associated with the domain of this operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.OperatorRangeMap | ( | self, | |
| args | |||
| ) |
OperatorRangeMap(PointRelaxation self) -> Map const Epetra_Map & Ifpack_PointRelaxation::OperatorRangeMap() const Returns the Epetra_Map object associated with the range of this operator.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.SetParameters | ( | self, | |
| args | |||
| ) |
SetParameters(PointRelaxation self, ParameterList List) -> int int Ifpack_PointRelaxation::SetParameters(Teuchos::ParameterList &List) Sets all the parameters for the preconditioner.
Reimplemented from PyTrilinos.IFPACK.Preconditioner.
| def PyTrilinos.IFPACK.PointRelaxation.SetUseTranspose | ( | self, | |
| args | |||
| ) |
SetUseTranspose(PointRelaxation self, bool UseTranspose_in) -> int virtual int Ifpack_PointRelaxation::SetUseTranspose(bool UseTranspose_in) This flag can be used to apply the preconditioner to the transpose of the input operator. Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.
Reimplemented from PyTrilinos.Epetra.Operator.
| def PyTrilinos.IFPACK.PointRelaxation.UseTranspose | ( | self, | |
| args | |||
| ) |
UseTranspose(PointRelaxation self) -> bool virtual bool Ifpack_PointRelaxation::UseTranspose() const Returns the current UseTranspose setting.
Reimplemented from PyTrilinos.Epetra.Operator.
1.7.6.1