43 #ifndef IFPACK_DENSECONTAINER_H 
   44 #define IFPACK_DENSECONTAINER_H 
   46 #include "Ifpack_ConfigDefs.h" 
   47 #include "Ifpack_Container.h" 
   48 #include "Epetra_SerialDenseMatrix.h" 
   49 #include "Epetra_SerialDenseSolver.h" 
   50 #include "Epetra_IntSerialDenseVector.h" 
  125     NumRows_(NumRows_in),
 
  126     NumVectors_(NumVectors_in),
 
  127     KeepNonFactoredMatrix_(false),
 
  128     IsInitialized_(false),
 
  132     ApplyInverseFlops_(0.0)
 
  144     if (KeepNonFactoredMatrix_)
 
  169     if (KeepNonFactoredMatrix_)
 
  194     if (NumVectors_ == NumVectors_in)
 
  197     NumVectors_ = NumVectors_in;
 
  198     IFPACK_CHK_ERR(RHS_.
Reshape(NumRows_,NumVectors_));
 
  199     IFPACK_CHK_ERR(LHS_.
Reshape(NumRows_,NumVectors_));
 
  201     for (
int i = 0 ; i < NumRows_ ; ++i)
 
  202       for (
int j = 0 ; j < NumVectors_ ; ++j) {
 
  208        IFPACK_CHK_ERR(Solver_.
SetVectors(LHS_,RHS_));
 
  214   virtual double& 
LHS(
const int i, 
const int Vector = 0);
 
  217   virtual double& 
RHS(
const int i, 
const int Vector = 0);
 
  229   virtual int& 
ID(
const int i);
 
  244     return(IsInitialized_);
 
  256     return(Label_.c_str());
 
  262     KeepNonFactoredMatrix_ = flag;
 
  269     return(KeepNonFactoredMatrix_);
 
  293     return(NonFactoredMatrix_);
 
  326     return(ComputeFlops_);
 
  336     return(ApplyInverseFlops_);
 
  340   virtual std::ostream& 
Print(std::ostream& os) 
const;
 
  364   bool KeepNonFactoredMatrix_;
 
  373   double ComputeFlops_;
 
  377   double ApplyInverseFlops_;
 
Ifpack_DenseContainer(const int NumRows_in, const int NumVectors_in=1)
Default constructor. 
virtual double ComputeFlops() const 
Returns the flops in Compute(). 
virtual const Epetra_IntSerialDenseVector & ID() const 
Returns the integer dense vector of IDs. 
virtual int Compute(const Epetra_RowMatrix &Matrix_in)
Finalizes the linear system matrix and prepares for the application of the inverse. 
virtual const char * Label() const 
Returns the label of this container. 
virtual std::ostream & Print(std::ostream &os) const 
Prints basic information on iostream. This function is used by operator<<. 
virtual int SetNumVectors(const int NumVectors_in)
Sets the number of vectors for LHS/RHS. 
virtual int SetMatrixElement(const int row, const int col, const double value)
Set the matrix element (row,col) to value. 
virtual int SetKeepNonFactoredMatrix(const bool flag)
If flag is true, keeps a copy of the non-factored matrix. 
virtual const Epetra_SerialDenseMatrix & Matrix() const 
Returns the dense matrix or its factors. 
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
virtual int NumVectors() const 
Returns the number of vectors in LHS/RHS. 
virtual double InitializeFlops() const 
Returns the flops in Initialize(). 
virtual double ApplyInverseFlops() const 
Returns the flops in ApplyInverse(). 
virtual bool IsInitialized() const 
Returns true is the container has been successfully initialized. 
virtual ~Ifpack_DenseContainer()
Destructor. 
Ifpack_DenseContainer(const Ifpack_DenseContainer &rhs)
Copy constructor. 
virtual int SetParameters(Teuchos::ParameterList &)
Sets all necessary parameters. 
virtual double & LHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of LHS. 
virtual int & ID(const int i)
Returns the ID associated to local row i. 
virtual int Apply()
Apply the matrix to RHS, results are stored in LHS. 
virtual bool IsComputed() const 
Returns true is the container has been successfully computed. 
virtual const Epetra_SerialDenseMatrix & LHS() const 
Returns the dense vector containing the LHS. 
Ifpack_DenseContainer & operator=(const Ifpack_DenseContainer &rhs)
Operator=. 
virtual int Initialize()
Initialize the container. 
int Reshape(int NumRows, int NumCols)
virtual int ApplyInverse()
Apply the inverse of the matrix to RHS, results are stored in LHS. 
virtual int NumRows() const 
Returns the number of rows of the matrix and LHS/RHS. 
Ifpack_DenseContainer: a class to define containers for dense matrices. 
Ifpack_Container: a pure virtual class for creating and solving local linear problems. 
virtual double & RHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of RHS. 
virtual double ApplyFlops() const 
Returns the flops in Apply(). 
virtual const Epetra_SerialDenseMatrix & NonFactoredMatrix() const 
Returns the non-factored dense matrix (only if stored). 
virtual bool KeepNonFactoredMatrix() const 
Returns KeepNonFactoredMatrix_. 
virtual const Epetra_SerialDenseMatrix & RHS() const 
Returns the dense vector containing the RHS.