46 #ifndef EPETRA_OSKIMATRIX_H
47 #define EPETRA_OSKIMATRIX_H
56 #include "oski/oski.h"
61 class Teuchos_ParameterList;
227 double Beta = 0.0)
const;
259 double Beta = 0.0)
const;
339 double Beta = 0.0)
const;
364 double Beta = 0.0)
const;
391 double Zeta = 0.0)
const;
418 double Zeta = 0.0)
const;
445 double Beta = 0.0)
const;
468 double Beta = 0.0)
const;
Epetra_MultiVector: A class for constructing and using dense multi-vectors, vectors and matrices in p...
int MatTransMatMultiply(bool ATA, const Epetra_Vector &x, Epetra_Vector &y, Epetra_Vector *t, double Alpha=1.0, double Beta=0.0) const
Performs two matrix vector multiplies of y = Alpha*this^TransA*this*x + Beta*y or y = Alpha*this*this...
int IsMatrixTransformed() const
Returns 1 if the matrix has been reordered by tuning, and 0 if it has not been.
Epetra_OskiPermutation: A class for storing the permutation performed on a Epetra_OskiMatrix.
Epetra_OskiMultiVector: A class for constructing and using dense Oski multi-vectors on a single proce...
const Epetra_OskiPermutation & ViewRowPermutation() const
Returns a read only row/left permutation of the Matrix.
int ExtractDiagonalCopy(Epetra_Vector &Diagonal) const
Returns a copy of the main diagonal in a user-provided vector.
int SetHintMatTransMatMultiply(bool ATA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, const Epetra_OskiMultiVector &Intermediate, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a two matrix-vector multiplies that are composed used by OskiTuneMat to ...
Epetra_OskiVector: A class for constructing and using dense OSKI vectors on a single processor or a s...
char * GetMatrixTransforms() const
Returns a string holding the transformations performed on the matrix when it was tuned.
Epetra_Vector: A class for constructing and using dense vectors on a parallel computer.
const Epetra_CrsMatrix * Epetra_View_
int ReplaceDiagonalValues(const Epetra_OskiVector &Diagonal)
Replaces diagonal values of the matrix with those in the user-provided vector.
Epetra_OskiMatrix: A class for constructing and using OSKI Matrices within Epetra. For information on known issues with OSKI see the detailed description.
int SetHintMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a matrix-vector multiply used by OskiTuneMat to optimize the data struct...
int MatPowMultiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y, Epetra_MultiVector &T, int Power=2, double Alpha=1.0, double Beta=0.0) const
Performs a matrix vector multiply of y = Alpha*(this^TransA)^Power*x + Beta*y. This is not implemente...
int TuneMatrix()
Tunes the matrix multiply if its deemed profitable.
int MultiplyAndMatTransMultiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y, const Epetra_Vector &w, Epetra_Vector &z, double Alpha=1.0, double Beta=0.0, double Omega=1.0, double Zeta=0.0) const
Performs the two matrix vector multiplies of y = Alpha*this*x + Beta*y and z = Omega*this^TransA*w + ...
Epetra_OskiMatrix(const Epetra_OskiMatrix &Source)
Copy constructor.
int SetHintMultiplyAndMatTransMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, double Omega, const Epetra_OskiMultiVector &InVec2, double Zeta, const Epetra_OskiMultiVector &OutVec2, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing two matrix-vector multiplies used by OskiTuneMat to optimize the data st...
int Solve(bool Upper, bool TransA, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const
Performs a triangular solve of y = (this^TransA)^-1*x where this is a triangular matrix.
int SetHint(const Teuchos::ParameterList &List)
Stores the hints in List in the matrix structure.
virtual ~Epetra_OskiMatrix()
Destructor.
int SetHintSolve(bool TransA, double Alpha, const Epetra_OskiMultiVector &Vector, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a triangular solve used by OskiTuneMat to optimize the data structure st...
int SumIntoMyValues(int MyRow, int NumEntries, double *Values, int *Indices)
Add this list of entries to existing values for a given local row of the matrix. WARNING: this could ...
const Epetra_OskiPermutation & ViewColumnPermutation() const
Returns a read only column/right permutation of the Matrix.
int SetHintPowMultiply(bool TransA, double Alpha, const Epetra_OskiMultiVector &InVec, double Beta, const Epetra_OskiMultiVector &OutVec, const Epetra_OskiMultiVector &Intermediate, int Power, int NumCalls, const Teuchos::ParameterList &List)
Workload hints for computing a matrix-vector multiply performed Power times used by OskiTuneMat to op...
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed r...
int Multiply(bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
Performs a matrix vector multiply of y = this^TransA*x.
int ReplaceMyValues(int MyRow, int NumEntries, double *Values, int *Indices)
Replace current values with this list of entries for a given local row of the matrix. Warning this could be expensive.
int ApplyMatrixTransforms(const char *Transforms)
Replaces the current data structure of the matrix with the one specified in Transforms.
const Epetra_OskiMatrix & ViewTransformedMat() const
Returns the transformed version of InMat if InMat has been transformed. If InMat has not been transfo...