Ifpack Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | List of all members
Ifpack_AdditiveSchwarz< T > Class Template Reference

Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix's. More...

#include <Ifpack_AdditiveSchwarz.h>

Inheritance diagram for Ifpack_AdditiveSchwarz< T >:
Inheritance graph
[legend]

Public Member Functions

virtual bool IsInitialized () const
 Returns true if the preconditioner has been successfully initialized. More...
 
virtual bool IsComputed () const
 Returns true if the preconditioner has been successfully computed. More...
 
virtual int SetParameters (Teuchos::ParameterList &List)
 Sets the parameters. More...
 
 Ifpack_AdditiveSchwarz (Epetra_RowMatrix *Matrix_in, int OverlapLevel_in=0)
 Ifpack_AdditiveSchwarz constructor with given Epetra_RowMatrix. More...
 
virtual ~Ifpack_AdditiveSchwarz ()
 Destructor. More...
 
virtual int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied (not implemented). More...
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Applies the matrix to X, returns the result in Y. More...
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Applies the preconditioner to X, returns the result in Y. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix (not implemented) More...
 
virtual const char * Label () const
 Returns a character string describing the operator. More...
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 
virtual bool HasNormInf () const
 Returns true if the this object can provide an approximate Inf-norm, false otherwise. More...
 
virtual const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
virtual const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
virtual const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator. More...
 
virtual int Initialize ()
 Initialized the preconditioner. More...
 
virtual int Compute ()
 Computes the preconditioner. More...
 
virtual double Condest (const Ifpack_CondestType CT=Ifpack_Cheap, const int MaxIters=1550, const double Tol=1e-9, Epetra_RowMatrix *Matrix_in=0)
 Computes the estimated condition number and returns its value. More...
 
virtual double Condest () const
 Returns the estimated condition number, or -1.0 if not computed. More...
 
virtual const Epetra_RowMatrixMatrix () const
 Returns a refernence to the internally stored matrix. More...
 
virtual bool IsOverlapping () const
 Returns true is an overlapping matrix is present. More...
 
virtual std::ostream & Print (std::ostream &) const
 Prints major information about this preconditioner. More...
 
virtual const T * Inverse () const
 
virtual int NumInitialize () const
 Returns the number of calls to Initialize(). More...
 
virtual int NumCompute () const
 Returns the number of calls to Compute(). More...
 
virtual int NumApplyInverse () const
 Returns the number of calls to ApplyInverse(). More...
 
virtual double InitializeTime () const
 Returns the time spent in Initialize(). More...
 
virtual double ComputeTime () const
 Returns the time spent in Compute(). More...
 
virtual double ApplyInverseTime () const
 Returns the time spent in ApplyInverse(). More...
 
virtual double InitializeFlops () const
 Returns the number of flops in the initialization phase. More...
 
virtual double ComputeFlops () const
 Returns the number of flops in the computation phase. More...
 
virtual double ApplyInverseFlops () const
 Returns the number of flops in the application of the preconditioner. More...
 
virtual int OverlapLevel () const
 Returns the level of overlap. More...
 
virtual const
Teuchos::ParameterList
List () const
 Returns a reference to the internally stored list. More...
 
 Ifpack_AdditiveSchwarz (const Ifpack_AdditiveSchwarz &RHS)
 Copy constructor (should never be used) More...
 
int Setup ()
 Sets up the localized matrix and the singleton filter. More...
 
Teuchos::RefCountPtr< const
Epetra_RowMatrix
Matrix_
 Pointers to the matrix to be preconditioned. More...
 
Teuchos::RefCountPtr
< Ifpack_OverlappingRowMatrix
OverlappingMatrix_
 Pointers to the overlapping matrix. More...
 
Teuchos::RefCountPtr
< Ifpack_LocalFilter
LocalizedMatrix_
 Localized version of Matrix_ or OverlappingMatrix_. More...
 
std::string Label_
 Contains the label of this object. More...
 
bool IsInitialized_
 If true, the preconditioner has been successfully initialized. More...
 
bool IsComputed_
 If true, the preconditioner has been successfully computed. More...
 
bool UseTranspose_
 If true, solve with the transpose (not supported by all solvers). More...
 
bool IsOverlapping_
 If true, overlapping is used. More...
 
int OverlapLevel_
 Level of overlap among the processors. More...
 
Teuchos::ParameterList List_
 Stores a copy of the list given in SetParameters() More...
 
Epetra_CombineMode CombineMode_
 Combine mode for off-process elements (only if overlap is used) More...
 
double Condest_
 Contains the estimated condition number. More...
 
bool ComputeCondest_
 If true, compute the condition number estimate each time Compute() is called. More...
 
bool UseReordering_
 If true, reorder the local matrix. More...
 
std::string ReorderingType_
 Type of reordering of the local matrix. More...
 
Teuchos::RefCountPtr
< Ifpack_Reordering
Reordering_
 Pointer to a reordering object. More...
 
Teuchos::RefCountPtr
< Ifpack_ReorderFilter
ReorderedLocalizedMatrix_
 Pointer to the reorderd matrix. More...
 
bool FilterSingletons_
 Filter for singletons. More...
 
Teuchos::RefCountPtr
< Ifpack_SingletonFilter
SingletonFilter_
 filtering object. More...
 
int NumInitialize_
 Contains the number of successful calls to Initialize(). More...
 
int NumCompute_
 Contains the number of successful call to Compute(). More...
 
int NumApplyInverse_
 Contains the number of successful call to ApplyInverse(). More...
 
double InitializeTime_
 Contains the time for all successful calls to Initialize(). More...
 
double ComputeTime_
 Contains the time for all successful calls to Compute(). More...
 
double ApplyInverseTime_
 Contains the time for all successful calls to ApplyInverse(). More...
 
double InitializeFlops_
 Contains the number of flops for Initialize(). More...
 
double ComputeFlops_
 Contains the number of flops for Compute(). More...
 
double ApplyInverseFlops_
 Contain sthe number of flops for ApplyInverse(). More...
 
Teuchos::RefCountPtr< Epetra_TimeTime_
 Object used for timing purposes. More...
 
Teuchos::RefCountPtr< T > Inverse_
 Pointer to the local solver. More...
 

Detailed Description

template<typename T>
class Ifpack_AdditiveSchwarz< T >

Ifpack_AdditiveSchwarz: a class to define Additive Schwarz preconditioners of Epetra_RowMatrix's.

Class Ifpack_AdditiveSchwarz enables the construction of Additive Schwarz (one-level overlapping domain decomposition) preconditioners, for a given Epetra_RowMatrix. Ifpack_AdditiveSchwarz is derived from Ifpack_Preconditioner, itself derived from Epetra_Operator. An application of the Additive Schwarz preconditioner can be obtained by calling method ApplyInverse().

One-level overlapping domain decomposition preconditioners use local solvers, of Dirichlet type. This means that the inverse of the local matrix (with minimal or wider overlap) is applied to the residual to be preconditioned.

The preconditioner can be written as:

\[ P_{AS}^{-1} = \sum_{i=1}^M P_i A_i^{-1} R_i , \]

where $M$ is the number of subdomains (that is, the number of processors in the computation), $R_i$ is an operator that restricts the global vector to the vector lying on subdomain $i$, $P_i$ is the prolongator operator, and

\[ A_i = R_i A P_i. \]

The construction of Schwarz preconditioners is mainly composed by two steps:

The local matrix $A_i$ can be filtered, to eliminate singletons, and reordered. At the present time, RCM and METIS can be used to reorder the local matrix.

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

Author
Marzio Sala, SNL 9214.
Date
Last modified on 22-Jan-05.

Definition at line 142 of file Ifpack_AdditiveSchwarz.h.

Constructor & Destructor Documentation

template<typename T >
Ifpack_AdditiveSchwarz< T >::Ifpack_AdditiveSchwarz ( Epetra_RowMatrix Matrix_in,
int  OverlapLevel_in = 0 
)

Ifpack_AdditiveSchwarz constructor with given Epetra_RowMatrix.

Creates an Ifpack_AdditiveSchwarz preconditioner with overlap. To use minimal-overlap, OverlappingMatrix is omitted (as defaulted to 0).

Parameters
Matrix- (In) Pointer to matrix to be preconditioned
OverlappingMatrix- (In) Pointer to the matrix extended with the desired level of overlap.

Definition at line 509 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual Ifpack_AdditiveSchwarz< T >::~Ifpack_AdditiveSchwarz ( )
inlinevirtual

Destructor.

Definition at line 163 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Ifpack_AdditiveSchwarz< T >::Ifpack_AdditiveSchwarz ( const Ifpack_AdditiveSchwarz< T > &  RHS)
inlineprotected

Copy constructor (should never be used)

Definition at line 382 of file Ifpack_AdditiveSchwarz.h.

Member Function Documentation

template<typename T >
int Ifpack_AdditiveSchwarz< T >::SetUseTranspose ( bool  UseTranspose_in)
virtual

If set true, transpose of this operator will be applied (not implemented).

This flag allows the transpose of the given operator to be used

implicitly.

Parameters
UseTranspose_in- (In) If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.

Implements Epetra_Operator.

Definition at line 931 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the matrix to X, returns the result in Y.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y-(Out) A Epetra_MultiVector of dimension NumVectors containing the result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 946 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Applies the preconditioner to X, returns the result in Y.

Parameters
X- (In) A Epetra_MultiVector of dimension NumVectors to be preconditioned.
Y-(Out) A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.
Warning
In order to work with AztecOO, any implementation of this method must support the case where X and Y are the same object.

Implements Ifpack_Preconditioner.

Definition at line 1004 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
double Ifpack_AdditiveSchwarz< T >::NormInf ( ) const
virtual

Returns the infinity norm of the global matrix (not implemented)

Implements Epetra_Operator.

Definition at line 954 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
const char * Ifpack_AdditiveSchwarz< T >::Label ( ) const
virtual

Returns a character string describing the operator.

Implements Epetra_Operator.

Definition at line 961 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
bool Ifpack_AdditiveSchwarz< T >::UseTranspose ( ) const
virtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 968 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
bool Ifpack_AdditiveSchwarz< T >::HasNormInf ( ) const
virtual

Returns true if the this object can provide an approximate Inf-norm, false otherwise.

Implements Epetra_Operator.

Definition at line 975 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
const Epetra_Comm & Ifpack_AdditiveSchwarz< T >::Comm ( ) const
virtual

Returns a pointer to the Epetra_Comm communicator associated with this operator.

Implements Epetra_Operator.

Definition at line 982 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
const Epetra_Map & Ifpack_AdditiveSchwarz< T >::OperatorDomainMap ( ) const
virtual

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

Implements Epetra_Operator.

Definition at line 989 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
const Epetra_Map & Ifpack_AdditiveSchwarz< T >::OperatorRangeMap ( ) const
virtual

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

Implements Epetra_Operator.

Definition at line 996 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual bool Ifpack_AdditiveSchwarz< T >::IsInitialized ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully initialized.

Implements Ifpack_Preconditioner.

Definition at line 236 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual bool Ifpack_AdditiveSchwarz< T >::IsComputed ( ) const
inlinevirtual

Returns true if the preconditioner has been successfully computed.

Implements Ifpack_Preconditioner.

Definition at line 242 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::SetParameters ( Teuchos::ParameterList List)
virtual

Sets the parameters.

Sets the parameter for the additive Schwarz preconditioner, as well as for all the preconditioners that may need to be defined on each subblock. Parameters accepted by List are:

  • "schwarz: combine mode" : It must be an Epetra_CombineMode. Default: Zero. It Can be assume of the following values:
    • Add: Components on the receiving processor will be added together;
    • Zero: Off-processor components will be ignored;
    • Insert: Off-processor components will be inserted into locations on receiving processor replacing existing values.
    • Average: Off-processor components will be averaged with existing;
    • AbsMax: Magnitudes of Off-processor components will be maxed with magnitudes of existing components on the receiving processor.
  • "schwarz: compute condest" : if true, Compute() will estimate the condition number of the preconditioner. Default: true.

Implements Ifpack_Preconditioner.

Definition at line 719 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::Initialize ( )
virtual

Initialized the preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 797 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::Compute ( )
virtual

Computes the preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 889 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
double Ifpack_AdditiveSchwarz< T >::Condest ( const Ifpack_CondestType  CT = Ifpack_Cheap,
const int  MaxIters = 1550,
const double  Tol = 1e-9,
Epetra_RowMatrix Matrix_in = 0 
)
virtual

Computes the estimated condition number and returns its value.

Implements Ifpack_Preconditioner.

Definition at line 1224 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::Condest ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 286 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual const Epetra_RowMatrix& Ifpack_AdditiveSchwarz< T >::Matrix ( ) const
inlinevirtual

Returns a refernence to the internally stored matrix.

Implements Ifpack_Preconditioner.

Definition at line 292 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual bool Ifpack_AdditiveSchwarz< T >::IsOverlapping ( ) const
inlinevirtual

Returns true is an overlapping matrix is present.

Definition at line 298 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
std::ostream & Ifpack_AdditiveSchwarz< T >::Print ( std::ostream &  os) const
virtual

Prints major information about this preconditioner.

Implements Ifpack_Preconditioner.

Definition at line 1143 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual const T* Ifpack_AdditiveSchwarz< T >::Inverse ( ) const
inlinevirtual

Definition at line 306 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual int Ifpack_AdditiveSchwarz< T >::NumInitialize ( ) const
inlinevirtual

Returns the number of calls to Initialize().

Implements Ifpack_Preconditioner.

Definition at line 312 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual int Ifpack_AdditiveSchwarz< T >::NumCompute ( ) const
inlinevirtual

Returns the number of calls to Compute().

Implements Ifpack_Preconditioner.

Definition at line 318 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual int Ifpack_AdditiveSchwarz< T >::NumApplyInverse ( ) const
inlinevirtual

Returns the number of calls to ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 324 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::InitializeTime ( ) const
inlinevirtual

Returns the time spent in Initialize().

Implements Ifpack_Preconditioner.

Definition at line 330 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::ComputeTime ( ) const
inlinevirtual

Returns the time spent in Compute().

Implements Ifpack_Preconditioner.

Definition at line 336 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::ApplyInverseTime ( ) const
inlinevirtual

Returns the time spent in ApplyInverse().

Implements Ifpack_Preconditioner.

Definition at line 342 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::InitializeFlops ( ) const
inlinevirtual

Returns the number of flops in the initialization phase.

Implements Ifpack_Preconditioner.

Definition at line 348 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::ComputeFlops ( ) const
inlinevirtual

Returns the number of flops in the computation phase.

Implements Ifpack_Preconditioner.

Definition at line 353 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual double Ifpack_AdditiveSchwarz< T >::ApplyInverseFlops ( ) const
inlinevirtual

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

Implements Ifpack_Preconditioner.

Definition at line 358 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual int Ifpack_AdditiveSchwarz< T >::OverlapLevel ( ) const
inlinevirtual

Returns the level of overlap.

Definition at line 364 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
virtual const Teuchos::ParameterList& Ifpack_AdditiveSchwarz< T >::List ( ) const
inlinevirtual

Returns a reference to the internally stored list.

Definition at line 370 of file Ifpack_AdditiveSchwarz.h.

template<typename T >
int Ifpack_AdditiveSchwarz< T >::Setup ( )
protected

Sets up the localized matrix and the singleton filter.

Definition at line 558 of file Ifpack_AdditiveSchwarz.h.

Member Data Documentation

template<typename T>
Teuchos::RefCountPtr<const Epetra_RowMatrix> Ifpack_AdditiveSchwarz< T >::Matrix_
protected

Pointers to the matrix to be preconditioned.

Definition at line 393 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Ifpack_OverlappingRowMatrix> Ifpack_AdditiveSchwarz< T >::OverlappingMatrix_
protected

Pointers to the overlapping matrix.

Definition at line 395 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Ifpack_LocalFilter> Ifpack_AdditiveSchwarz< T >::LocalizedMatrix_
protected

Localized version of Matrix_ or OverlappingMatrix_.

Definition at line 424 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
std::string Ifpack_AdditiveSchwarz< T >::Label_
protected

Contains the label of this object.

Definition at line 441 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::IsInitialized_
protected

If true, the preconditioner has been successfully initialized.

Definition at line 443 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::IsComputed_
protected

If true, the preconditioner has been successfully computed.

Definition at line 445 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::UseTranspose_
protected

If true, solve with the transpose (not supported by all solvers).

Definition at line 447 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::IsOverlapping_
protected

If true, overlapping is used.

Definition at line 449 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
int Ifpack_AdditiveSchwarz< T >::OverlapLevel_
protected

Level of overlap among the processors.

Definition at line 451 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::ParameterList Ifpack_AdditiveSchwarz< T >::List_
protected

Stores a copy of the list given in SetParameters()

Definition at line 453 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Epetra_CombineMode Ifpack_AdditiveSchwarz< T >::CombineMode_
protected

Combine mode for off-process elements (only if overlap is used)

Definition at line 455 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::Condest_
protected

Contains the estimated condition number.

Definition at line 457 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::ComputeCondest_
protected

If true, compute the condition number estimate each time Compute() is called.

Definition at line 459 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::UseReordering_
protected

If true, reorder the local matrix.

Definition at line 461 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
std::string Ifpack_AdditiveSchwarz< T >::ReorderingType_
protected

Type of reordering of the local matrix.

Definition at line 463 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Ifpack_Reordering> Ifpack_AdditiveSchwarz< T >::Reordering_
protected

Pointer to a reordering object.

Definition at line 465 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Ifpack_ReorderFilter> Ifpack_AdditiveSchwarz< T >::ReorderedLocalizedMatrix_
protected

Pointer to the reorderd matrix.

Definition at line 467 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
bool Ifpack_AdditiveSchwarz< T >::FilterSingletons_
protected

Filter for singletons.

Definition at line 469 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Ifpack_SingletonFilter> Ifpack_AdditiveSchwarz< T >::SingletonFilter_
protected

filtering object.

Definition at line 471 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
int Ifpack_AdditiveSchwarz< T >::NumInitialize_
protected

Contains the number of successful calls to Initialize().

Definition at line 473 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
int Ifpack_AdditiveSchwarz< T >::NumCompute_
protected

Contains the number of successful call to Compute().

Definition at line 475 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
int Ifpack_AdditiveSchwarz< T >::NumApplyInverse_
mutableprotected

Contains the number of successful call to ApplyInverse().

Definition at line 477 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::InitializeTime_
protected

Contains the time for all successful calls to Initialize().

Definition at line 479 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::ComputeTime_
protected

Contains the time for all successful calls to Compute().

Definition at line 481 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::ApplyInverseTime_
mutableprotected

Contains the time for all successful calls to ApplyInverse().

Definition at line 483 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::InitializeFlops_
protected

Contains the number of flops for Initialize().

Definition at line 485 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::ComputeFlops_
protected

Contains the number of flops for Compute().

Definition at line 487 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
double Ifpack_AdditiveSchwarz< T >::ApplyInverseFlops_
mutableprotected

Contain sthe number of flops for ApplyInverse().

Definition at line 489 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<Epetra_Time> Ifpack_AdditiveSchwarz< T >::Time_
protected

Object used for timing purposes.

Definition at line 491 of file Ifpack_AdditiveSchwarz.h.

template<typename T>
Teuchos::RefCountPtr<T> Ifpack_AdditiveSchwarz< T >::Inverse_
protected

Pointer to the local solver.

Definition at line 493 of file Ifpack_AdditiveSchwarz.h.


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