43 #ifndef IFPACK_TRIDICONTAINER_H
44 #define IFPACK_TRIDICONTAINER_H
50 #include "Epetra_IntSerialDenseVector.h"
51 #include "Epetra_SerialDenseVector.h"
202 for (
int i = 0 ; i <
NumRows_ ; ++i)
215 virtual double&
LHS(
const int i,
const int Vector = 0);
218 virtual double&
RHS(
const int i,
const int Vector = 0);
230 virtual int&
ID(
const int i);
341 virtual std::ostream&
Print(std::ostream& os)
const;
virtual int NumRows() const
Returns the number of rows of the matrix and LHS/RHS.
virtual std::ostream & Print(std::ostream &os) const
Prints basic information on iostream. This function is used by operator<<.
std::string Label_
Label for this object.
Ifpack_TriDiContainer(const Ifpack_TriDiContainer &rhs)
Copy constructor.
double ApplyInverseFlops_
Flops in ApplyInverse().
virtual int NumVectors() const
Returns the number of vectors in LHS/RHS.
Ifpack_TriDiContainer: a class to define containers for dense matrices.
virtual const Epetra_IntSerialDenseVector & ID() const
Returns the integer dense vector of IDs.
virtual int SetMatrixElement(const int row, const int col, const double value)
Set the matrix element (row,col) to value.
bool IsInitialized_
If true, the container has been successfully initialized.
bool KeepNonFactoredMatrix_
If true, keeps a copy of the non-factored matrix.
virtual int Extract(const Epetra_RowMatrix &Matrix_in)
Extract the submatrices identified by the ID set int ID().
virtual const Epetra_SerialDenseMatrix & RHS() const
Returns the dense vector containing the RHS.
virtual double InitializeFlops() const
Returns the flops in Initialize().
Epetra_SerialDenseMatrix RHS_
SerialDense vector representing the RHS.
virtual double & LHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of LHS.
Ifpack_SerialTriDiSolver: A class for solving TriDi linear problems.
virtual const char * Label() const
Returns the label of this container.
virtual const Epetra_SerialDenseMatrix & LHS() const
Returns the dense vector containing the LHS.
int NumVectors_
Number of vectors in the container.
virtual int SetParameters(Teuchos::ParameterList &)
Sets all necessary parameters.
virtual int SetKeepNonFactoredMatrix(const bool flag)
If flag is true, keeps a copy of the non-factored matrix.
Ifpack_SerialTriDiSolver Solver_
TriDi solver (solution will be get using LAPACK).
int NumRows_
Number of rows in the container.
virtual int & ID(const int i)
Returns the ID associated to local row i.
virtual const Ifpack_SerialTriDiMatrix & Matrix() const
Returns the dense matrix or its factors.
virtual double & RHS(const int i, const int Vector=0)
Returns the i-th component of the vector Vector of RHS.
virtual ~Ifpack_TriDiContainer()
Destructor.
virtual double ComputeFlops() const
Returns the flops in Compute().
Epetra_IntSerialDenseVector ID_
Sets of local rows.
virtual double ApplyFlops() const
Returns the flops in Apply().
virtual int SetNumVectors(const int NumVectors_in)
Sets the number of vectors for LHS/RHS.
virtual const Ifpack_SerialTriDiMatrix & NonFactoredMatrix() const
Returns the non-factored dense matrix (only if stored).
Ifpack_TriDiContainer & operator=(const Ifpack_TriDiContainer &rhs)
Operator=.
int Reshape(int NumRows, int NumCols)
virtual int Apply()
Apply the matrix to RHS, results are stored in LHS.
double ApplyFlops_
Flops in Apply().
virtual double ApplyInverseFlops() const
Returns the flops in ApplyInverse().
virtual int Compute(const Epetra_RowMatrix &Matrix_in)
Finalizes the linear system matrix and prepares for the application of the inverse.
Ifpack_SerialTriDiMatrix: A class for constructing and using real double precision general TriDi matr...
virtual int ApplyInverse()
Apply the inverse of the matrix to RHS, results are stored in LHS.
Ifpack_SerialTriDiMatrix Matrix_
TriDi matrix.
virtual int Initialize()
Initialize the container.
Ifpack_SerialTriDiMatrix NonFactoredMatrix_
TriDi matrix, that contains the non-factored matrix.
int SetVectors(Epetra_SerialDenseMatrix &X, Epetra_SerialDenseMatrix &B)
Sets the pointers for left and right hand side vector(s).
Ifpack_Container: a pure virtual class for creating and solving local linear problems.
bool IsComputed_
If true, the container has been successfully computed.
#define IFPACK_CHK_ERR(ifpack_err)
double ComputeFlops_
Flops in Compute().
virtual bool IsInitialized() const
Returns true is the container has been successfully initialized.
virtual bool KeepNonFactoredMatrix() const
Returns KeepNonFactoredMatrix_.
virtual bool IsComputed() const
Returns true is the container has been successfully computed.
Epetra_SerialDenseMatrix LHS_
SerialDense vector representing the LHS.
Ifpack_TriDiContainer(const int NumRows_in, const int NumVectors_in=1)
Default constructor.