43 #ifndef IFPACK_CONTAINER_H 
   44 #define IFPACK_CONTAINER_H 
  106   virtual int NumRows() 
const = 0;
 
  115   virtual double& 
LHS(
const int i, 
const int Vector = 0) = 0;
 
  118   virtual double& 
RHS(
const int i, 
const int Vector = 0) = 0;
 
  130   virtual int& 
ID(
const int i) = 0;
 
  134                                const double value) = 0;
 
  152   virtual int Apply() = 0;
 
  158   virtual const char* 
Label() 
const = 0;
 
  173   virtual std::ostream& 
Print(std::ostream& os) 
const = 0;
 
  178   return(obj.
Print(os));
 
  181 #endif // IFPACK_CONTAINER_H 
virtual int SetNumVectors(const int i)=0
Sets the number of vectors for LHS/RHS. 
 
virtual double InitializeFlops() const =0
Returns the flops in Initialize(). 
 
virtual int SetParameters(Teuchos::ParameterList &List)=0
Sets all necessary parameters. 
 
virtual ~Ifpack_Container()
Destructor. 
 
virtual bool IsInitialized() const =0
Returns true is the container has been successfully initialized. 
 
virtual double ApplyInverseFlops() const =0
Returns the flops in ApplyInverse(). 
 
virtual const char * Label() const =0
Returns the label of this container. 
 
virtual int Initialize()=0
Initializes the container, by performing all operations that only require matrix structure. 
 
virtual int Apply()=0
Apply the matrix to RHS, results are stored in LHS. 
 
virtual int SetMatrixElement(const int row, const int col, const double value)=0
Set the matrix element (row,col) to value. 
 
virtual int Compute(const Epetra_RowMatrix &A)=0
Finalizes the linear system matrix and prepares for the application of the inverse. 
 
Ifpack_Partitioner: A class to decompose local Ifpack_Graph's. 
 
virtual int ApplyInverse()=0
Apply the inverse of the matrix to RHS, results are stored in LHS. 
 
virtual double & RHS(const int i, const int Vector=0)=0
Returns the i-th component of the vector Vector of RHS. 
 
virtual double ApplyFlops() const =0
Returns the flops in Apply(). 
 
virtual std::ostream & Print(std::ostream &os) const =0
Prints out basic information about the container. 
 
virtual double ComputeFlops() const =0
Returns the flops in Compute(). 
 
std::ostream & operator<<(std::ostream &os, const Ifpack_Container &obj)
 
virtual int NumRows() const =0
Returns the number of rows of the matrix and LHS/RHS. 
 
Ifpack_Container: a pure virtual class for creating and solving local linear problems. 
 
virtual int NumVectors() const =0
Returns the number of vectors in LHS/RHS. 
 
virtual bool IsComputed() const =0
Returns true is the container has been successfully computed. 
 
virtual double & LHS(const int i, const int Vector=0)=0
Returns the i-th component of the vector Vector of LHS. 
 
virtual int & ID(const int i)=0
Returns the ID associated to local row i.