Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
List of all members
Epetra_OskiMatrix Class Reference

Epetra_OskiMatrix: A class for constructing and using OSKI Matrices within Epetra. For information on known issues with OSKI see the detailed description. More...

#include <Epetra_OskiMatrix.h>

Inheritance diagram for Epetra_OskiMatrix:
Inheritance graph
[legend]
Collaboration diagram for Epetra_OskiMatrix:
Collaboration graph
[legend]

Public Member Functions

Constructors/Destructor
 Epetra_OskiMatrix (const Epetra_OskiMatrix &Source)
 Copy constructor.
 
 Epetra_OskiMatrix (const Epetra_CrsMatrix &Source, const Teuchos::ParameterList &List)
 Constructor creates an Epetra_OskiMatrix from an Epetra_CrsMatrix. More...
 
virtual ~Epetra_OskiMatrix ()
 Destructor.
 
Extract/Replace Values
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. More...
 
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 be expensive. More...
 
int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 Returns a copy of the main diagonal in a user-provided vector. More...
 
int ReplaceDiagonalValues (const Epetra_OskiVector &Diagonal)
 Replaces diagonal values of the matrix with those in the user-provided vector. More...
 
Computational methods
int Multiply (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
 Performs a matrix vector multiply of y = this^TransA*x. More...
 
int Multiply (bool TransA, const Epetra_Vector &x, Epetra_Vector &y, double Alpha, double Beta=0.0) const
 Performs a matrix vector multiply of y = Alpha*this^TransA*x + Beta*y. More...
 
int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Performs a matrix multi-vector multiply of Y = this^TransA*X. More...
 
int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y, double Alpha, double Beta=0.0) const
 Performs a matrix multi-vector multiply of Y = Alpha*this^TransA*X + Beta*Y. More...
 
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. More...
 
int Solve (bool TransA, const Epetra_Vector &x, Epetra_Vector &y, double Alpha=1.0) const
 Performs a triangular solve of y = Alpha*(this^TransA)^-1*x where this is a triangular matrix. More...
 
int Solve (bool Upper, bool TransA, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Performs a triangular solve of Y = (this^TransA)^-1*X where this is a triangular matrix. More...
 
int Solve (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y, double Alpha=1.0) const
 Performs a triangular solve of Y = Alpha*(this^TransA)^-1*X where this is a triangular matrix. More...
 
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^TransA*x + Beta*y. More...
 
int MatTransMatMultiply (bool ATA, const Epetra_MultiVector &X, Epetra_MultiVector &Y, Epetra_MultiVector *T, double Alpha=1.0, double Beta=0.0) const
 Performs two matrix multi-vector multiplies of Y = Alpha*this^TransA*this*X + Beta*Y or Y = Alpha*this*this^TransA*X + Beta*Y. More...
 
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 + Zeta*z. More...
 
int MultiplyAndMatTransMultiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y, const Epetra_MultiVector &W, Epetra_MultiVector &Z, double Alpha=1.0, double Beta=0.0, double Omega=1.0, double Zeta=0.0) const
 Performs the two matrix multi-vector multiplies of Y = Alpha*this*X + Beta*Y and Z = Omega*this^TransA*W + Zeta*Z. More...
 
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 implemented as described in the detailed description. More...
 
int MatPowMultiply (bool TransA, const Epetra_Vector &x, Epetra_Vector &y, 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 implemented as described in the detailed description. More...
 
Tuning
int SetHint (const Teuchos::ParameterList &List)
 Stores the hints in List in the matrix structure. More...
 
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 structure storage, and the routine to compute the calculation. More...
 
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 storage, and the routine to compute the calculation. More...
 
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 optimize the data structure storage, and the routine to compute the calculation. More...
 
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 structure storage, and the routine to compute the calculation. More...
 
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 optimize the data structure storage and the routine to compute the calculation. More...
 
int TuneMatrix ()
 Tunes the matrix multiply if its deemed profitable. More...
 
Data Structure Transformation Methods
int IsMatrixTransformed () const
 Returns 1 if the matrix has been reordered by tuning, and 0 if it has not been.
 
const Epetra_OskiMatrixViewTransformedMat () const
 Returns the transformed version of InMat if InMat has been transformed. If InMat has not been transformed then the return will equal InMat.
 
const Epetra_OskiPermutationViewRowPermutation () const
 Returns a read only row/left permutation of the Matrix.
 
const Epetra_OskiPermutationViewColumnPermutation () const
 Returns a read only column/right permutation of the Matrix.
 
char * GetMatrixTransforms () const
 Returns a string holding the transformations performed on the matrix when it was tuned. More...
 
int ApplyMatrixTransforms (const char *Transforms)
 Replaces the current data structure of the matrix with the one specified in Transforms. More...
 
- Public Member Functions inherited from Epetra_CrsMatrix
void FusedImport (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 
void FusedExport (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 
void FusedImport (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Import *DomainImporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 
void FusedExport (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Export *DomainExporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const int *NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constructor with variable number of indices per row. More...
 
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constructor with fixed number of indices per row. More...
 
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, const int *NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constructor with variable number of indices per row. More...
 
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, int NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constuctor with fixed number of indices per row. More...
 
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_CrsGraph &Graph)
 Construct a matrix using an existing Epetra_CrsGraph object. More...
 
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false)
 Epetra CrsMatrix constructor that also fuses Import and FillComplete(). More...
 
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Import *DomainImporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 Epetra CrsMatrix constructor that also fuses Import and FillComplete(). More...
 
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false)
 Epetra CrsMatrix constructor that also fuses Ex[prt and FillComplete(). More...
 
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Export *DomainExporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator)
 Epetra CrsMatrix constructor that also fuses Ex[prt and FillComplete(). More...
 
 Epetra_CrsMatrix (const Epetra_CrsMatrix &Matrix)
 Copy constructor.
 
virtual ~Epetra_CrsMatrix ()
 Epetra_CrsMatrix Destructor.
 
Epetra_CrsMatrixoperator= (const Epetra_CrsMatrix &src)
 Assignment operator.
 
int PutScalar (double ScalarConstant)
 Initialize all values in the matrix with constant value. More...
 
int Scale (double ScalarConstant)
 Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A). More...
 
virtual int InsertGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices)
 Insert a list of elements in a given global row of the matrix. More...
 
virtual int InsertGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices)
 
virtual int InsertGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices)
 Insert a list of elements in a given global row of the matrix. More...
 
virtual int InsertGlobalValues (long long GlobalRow, int NumEntries, double *Values, long long *Indices)
 
virtual int ReplaceGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices)
 Replace specified existing values with this list of entries for a given global row of the matrix. More...
 
virtual int ReplaceGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices)
 
virtual int SumIntoGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices)
 Add this list of entries to existing values for a given global row of the matrix. More...
 
virtual int SumIntoGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices)
 
int InsertMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices)
 Insert a list of elements in a given local row of the matrix. More...
 
int InsertMyValues (int MyRow, int NumEntries, double *Values, int *Indices)
 Insert a list of elements in a given local row of the matrix. More...
 
int ReplaceMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices)
 Replace current values with this list of entries for a given local row of the matrix. More...
 
int SumIntoMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices)
 Add this list of entries to existing values for a given local row of the matrix. More...
 
int ReplaceDiagonalValues (const Epetra_Vector &Diagonal)
 Replaces diagonal values of the matrix with those in the user-provided vector. More...
 
int FillComplete (bool OptimizeDataStorage=true)
 Signal that data entry is complete. Perform transformations to local index space.
 
int FillComplete (const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, bool OptimizeDataStorage=true)
 Signal that data entry is complete. Perform transformations to local index space.
 
int OptimizeStorage ()
 Make consecutive row index sections contiguous, minimize internal storage used for constructing graph. More...
 
int MakeDataContiguous ()
 Eliminates memory that is used for construction. Make consecutive row index sections contiguous.
 
int ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified global row in user-provided arrays. More...
 
int ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values, long long *Indices) const
 
int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified local row in user-provided arrays. More...
 
int ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values) const
 Returns a copy of the specified global row values in user-provided array. More...
 
int ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values) const
 
int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values) const
 Returns a copy of the specified local row values in user-provided array. More...
 
int ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values, int *&Indices) const
 Returns a view of the specified global row values via pointers to internal data. More...
 
int ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values, long long *&Indices) const
 
int ExtractMyRowView (int MyRow, int &NumEntries, double *&Values, int *&Indices) const
 Returns a view of the specified local row values via pointers to internal data. More...
 
int ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values) const
 Returns a view of the specified global row values via pointers to internal data. More...
 
int ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values) const
 
int ExtractMyRowView (int MyRow, int &NumEntries, double *&Values) const
 Returns a view of the specified local row values via pointers to internal data. More...
 
int Multiply (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
 Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_Vector x in y. More...
 
int Multiply1 (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
 
int Multiply1 (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 
int Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const
 Returns the result of a local solve using the Epetra_CrsMatrix on a Epetra_Vector x in y. More...
 
int InvRowSums (Epetra_Vector &x) const
 Computes the inverse of the sum of absolute values of the rows of the Epetra_CrsMatrix, results returned in x. More...
 
int InvRowMaxs (Epetra_Vector &x) const
 Computes the inverse of the max of absolute values of the rows of the Epetra_CrsMatrix, results returned in x. More...
 
int LeftScale (const Epetra_Vector &x)
 Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x. More...
 
int InvColSums (Epetra_Vector &x) const
 Computes the inverse of the sum of absolute values of the columns of the Epetra_CrsMatrix, results returned in x. More...
 
int InvColMaxs (Epetra_Vector &x) const
 Computes the max of absolute values of the columns of the Epetra_CrsMatrix, results returned in x. More...
 
int RightScale (const Epetra_Vector &x)
 Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x. More...
 
bool Filled () const
 If FillComplete() has been called, this query returns true, otherwise it returns false.
 
bool StorageOptimized () const
 If OptimizeStorage() has been called, this query returns true, otherwise it returns false.
 
bool IndicesAreGlobal () const
 If matrix indices has not been transformed to local, this query returns true, otherwise it returns false.
 
bool IndicesAreLocal () const
 If matrix indices has been transformed to local, this query returns true, otherwise it returns false.
 
bool IndicesAreContiguous () const
 If matrix indices are packed into single array (done in OptimizeStorage()) return true, otherwise false.
 
bool LowerTriangular () const
 If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.
 
bool UpperTriangular () const
 If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.
 
bool NoDiagonal () const
 If matrix has no diagonal entries in global index space, this query returns true, otherwise it returns false.
 
double NormInf () const
 Returns the infinity norm of the global matrix.
 
double NormOne () const
 Returns the one norm of the global matrix.
 
double NormFrobenius () const
 Returns the frobenius norm of the global matrix.
 
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global matrix.
 
long long NumGlobalNonzeros64 () const
 
int NumGlobalRows () const
 Returns the number of global matrix rows.
 
long long NumGlobalRows64 () const
 
int NumGlobalCols () const
 Returns the number of global matrix columns.
 
long long NumGlobalCols64 () const
 
int NumGlobalDiagonals () const
 Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.
 
long long NumGlobalDiagonals64 () const
 
int NumMyNonzeros () const
 Returns the number of nonzero entries in the calling processor's portion of the matrix.
 
int NumMyRows () const
 Returns the number of matrix rows owned by the calling processor.
 
int NumMyCols () const
 Returns the number of entries in the set of column-indices that appear on this processor. More...
 
int NumMyDiagonals () const
 Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. More...
 
int NumGlobalEntries (long long Row) const
 Returns the current number of nonzero entries in specified global row on this processor.
 
int NumAllocatedGlobalEntries (int Row) const
 Returns the allocated number of nonzero entries in specified global row on this processor.
 
int MaxNumEntries () const
 Returns the maximum number of nonzero entries across all rows on this processor. More...
 
int GlobalMaxNumEntries () const
 Returns the maximum number of nonzero entries across all rows on all processors. More...
 
int NumMyEntries (int Row) const
 Returns the current number of nonzero entries in specified local row on this processor.
 
int NumAllocatedMyEntries (int Row) const
 Returns the allocated number of nonzero entries in specified local row on this processor.
 
int IndexBase () const
 Returns the index base for row and column indices for this graph. More...
 
long long IndexBase64 () const
 
bool StaticGraph ()
 Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable.
 
const Epetra_CrsGraphGraph () const
 Returns a reference to the Epetra_CrsGraph object associated with this matrix.
 
const Epetra_MapRowMap () const
 Returns the Epetra_Map object associated with the rows of this matrix.
 
int ReplaceRowMap (const Epetra_BlockMap &newmap)
 Replaces the current RowMap with the user-specified map object. More...
 
bool HaveColMap () const
 Returns true if we have a well-defined ColMap, and returns false otherwise. More...
 
int ReplaceColMap (const Epetra_BlockMap &newmap)
 Replaces the current ColMap with the user-specified map object. More...
 
int ReplaceDomainMapAndImporter (const Epetra_Map &NewDomainMap, const Epetra_Import *NewImporter)
 Replaces the current DomainMap & Importer with the user-specified map object. More...
 
int RemoveEmptyProcessesInPlace (const Epetra_BlockMap *NewMap)
 Remove processes owning zero rows from the Maps and their communicator. More...
 
const Epetra_MapColMap () const
 Returns the Epetra_Map object that describes the set of column-indices that appear in each processor's locally owned matrix rows. More...
 
const Epetra_MapDomainMap () const
 Returns the Epetra_Map object associated with the domain of this matrix operator. More...
 
const Epetra_MapRangeMap () const
 Returns the Epetra_Map object associated with the range of this matrix operator. More...
 
const Epetra_ImportImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations.
 
const Epetra_ExportExporter () const
 Returns the Epetra_Export object that contains the export operations for distributed operations.
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix.
 
int LRID (int GRID_in) const
 Returns the local row index for given global row index, returns -1 if no local row for this global row.
 
int LRID (long long GRID_in) const
 
int GRID (int LRID_in) const
 Returns the global row index for give local row index, returns IndexBase-1 if we don't have this local row.
 
long long GRID64 (int LRID_in) const
 
int LCID (int GCID_in) const
 Returns the local column index for given global column index, returns -1 if no local column for this global column. More...
 
int LCID (long long GCID_in) const
 
int GCID (int LCID_in) const
 Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column. More...
 
long long GCID64 (int LCID_in) const
 
bool MyGRID (int GRID_in) const
 Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false.
 
bool MyGRID (long long GRID_in) const
 
bool MyLRID (int LRID_in) const
 Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
 
bool MyGCID (int GCID_in) const
 Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false. More...
 
bool MyGCID (long long GCID_in) const
 
bool MyLCID (int LCID_in) const
 Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false. More...
 
bool MyGlobalRow (int GID) const
 Returns true of GID is owned by the calling processor, otherwise it returns false.
 
bool MyGlobalRow (long long GID) const
 
virtual void Print (std::ostream &os) const
 Print method.
 
const char * Label () const
 Returns a character string describing the operator.
 
int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied. More...
 
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
bool HasNormInf () const
 Returns true because this class can compute an Inf-norm.
 
bool UseTranspose () const
 Returns the current UseTranspose setting.
 
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this matrix operator.
 
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this matrix operator.
 
int NumMyRowEntries (int MyRow, int &NumEntries) const
 Return the current number of values stored for the specified local row. More...
 
const Epetra_BlockMapMap () const
 Map() method inherited from Epetra_DistObject.
 
const Epetra_MapRowMatrixRowMap () const
 Returns the Epetra_Map object associated with the rows of this matrix.
 
const Epetra_MapRowMatrixColMap () const
 Returns the Epetra_Map object associated with columns of this matrix.
 
const Epetra_ImportRowMatrixImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations.
 
double * operator[] (int Loc)
 Inlined bracket operator for fast access to data. (Const and Non-const versions) More...
 
double * operator[] (int Loc) const
 
int ExtractCrsDataPointers (int *&IndexOffset, int *&Indices, double *&Values_in) const
 Returns internal data pointers associated with Crs matrix format. More...
 
Epetra_IntSerialDenseVectorExpertExtractIndexOffset ()
 Returns a reference to the Epetra_IntSerialDenseVector used to hold the local IndexOffsets (CRS rowptr) More...
 
Epetra_IntSerialDenseVectorExpertExtractIndices ()
 Returns a reference to the Epetra_IntSerialDenseVector used to hold the local All_Indices (CRS colind) More...
 
double *& ExpertExtractValues ()
 Returns a reference to the double* used to hold the values array. More...
 
int ExpertStaticFillComplete (const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, const Epetra_Import *Importer=0, const Epetra_Export *Exporter=0, int NumMyDiagonals=-1)
 Performs a FillComplete on an object that aready has filled CRS data. More...
 
int ExpertMakeUniqueCrsGraphData ()
 Makes sure this matrix has a unique CrsGraphData object. More...
 
int SortGhostsAssociatedWithEachProcessor (bool Flag)
 Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order. More...
 
const Epetra_MapImportMap () const
 Use ColMap() instead.
 
int TransformToLocal ()
 Use FillComplete() instead.
 
int TransformToLocal (const Epetra_Map *DomainMap, const Epetra_Map *RangeMap)
 Use FillComplete(const Epetra_Map& DomainMap, const Epetra_Map& RangeMap) instead.
 
- Public Member Functions inherited from Epetra_DistObject
 Epetra_DistObject (const Epetra_BlockMap &Map)
 Basic Epetra_DistObject constuctor. More...
 
 Epetra_DistObject (const Epetra_BlockMap &Map, const char *const Label)
 
 Epetra_DistObject (const Epetra_DistObject &Source)
 Epetra_DistObject copy constructor.
 
virtual ~Epetra_DistObject ()
 Epetra_DistObject destructor.
 
int Import (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Imports an Epetra_DistObject using the Epetra_Import object. More...
 
int Import (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Imports an Epetra_DistObject using the Epetra_Export object. More...
 
int Export (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Exports an Epetra_DistObject using the Epetra_Import object. More...
 
int Export (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Exports an Epetra_DistObject using the Epetra_Export object. More...
 
const Epetra_BlockMapMap () const
 Returns the address of the Epetra_BlockMap for this multi-vector.
 
const Epetra_CommComm () const
 Returns the address of the Epetra_Comm for this multi-vector.
 
bool DistributedGlobal () const
 Returns true if this multi-vector is distributed global, i.e., not local replicated.
 
- Public Member Functions inherited from Epetra_Object
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const Epetra_Object &Object)
 Epetra_Object Copy Constructor. More...
 
virtual ~Epetra_Object ()
 Epetra_Object Destructor. More...
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method.
 
virtual void SetLabel (const char *const Label)
 Epetra_Object Label definition using char *. More...
 
- Public Member Functions inherited from Epetra_SrcDistObject
virtual ~Epetra_SrcDistObject ()
 Epetra_SrcDistObject destructor.
 
- Public Member Functions inherited from Epetra_CompObject
Epetra_CompObjectoperator= (const Epetra_CompObject &src)
 
 Epetra_CompObject ()
 Basic Epetra_CompObject constuctor.
 
 Epetra_CompObject (const Epetra_CompObject &Source)
 Epetra_CompObject copy constructor.
 
virtual ~Epetra_CompObject ()
 Epetra_CompObject destructor.
 
void SetFlopCounter (const Epetra_Flops &FlopCounter_in)
 Set the internal Epetra_Flops() pointer.
 
void SetFlopCounter (const Epetra_CompObject &CompObject)
 Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject.
 
void UnsetFlopCounter ()
 Set the internal Epetra_Flops() pointer to 0 (no flops counted).
 
Epetra_FlopsGetFlopCounter () const
 Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none.
 
void ResetFlops () const
 Resets the number of floating point operations to zero for this multi-vector.
 
double Flops () const
 Returns the number of floating point operations with this multi-vector.
 
void UpdateFlops (int Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (long int Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (long long Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (double Flops_in) const
 Increment Flop count for this object.
 
void UpdateFlops (float Flops_in) const
 Increment Flop count for this object.
 
- Public Member Functions inherited from Epetra_BLAS
 Epetra_BLAS (void)
 Epetra_BLAS Constructor. More...
 
 Epetra_BLAS (const Epetra_BLAS &BLAS)
 Epetra_BLAS Copy Constructor. More...
 
virtual ~Epetra_BLAS (void)
 Epetra_BLAS Destructor.
 
float ASUM (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS one norm function (SASUM).
 
double ASUM (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS one norm function (DASUM).
 
float DOT (const int N, const float *X, const float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (SDOT).
 
double DOT (const int N, const double *X, const double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS dot product function (DDOT).
 
float NRM2 (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS norm function (SNRM2).
 
double NRM2 (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS norm function (DNRM2).
 
void SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const
 Epetra_BLAS vector scale function (SSCAL)
 
void SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const
 Epetra_BLAS vector scale function (DSCAL)
 
void COPY (const int N, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector copy function (SCOPY)
 
void COPY (const int N, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector scale function (DCOPY)
 
int IAMAX (const int N, const float *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (ISAMAX)
 
int IAMAX (const int N, const double *X, const int INCX=1) const
 Epetra_BLAS arg maximum of absolute value function (IDAMAX)
 
void AXPY (const int N, const float ALPHA, const float *X, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (SAXPY)
 
void AXPY (const int N, const double ALPHA, const double *X, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS vector update function (DAXPY)
 
void GEMV (const char TRANS, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *X, const float BETA, float *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (SGEMV)
 
void GEMV (const char TRANS, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *X, const double BETA, double *Y, const int INCX=1, const int INCY=1) const
 Epetra_BLAS matrix-vector multiply function (DGEMV)
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (SGEMM)
 
void GEMM (const char TRANSA, const char TRANSB, const int M, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DGEMM)
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const float ALPHA, const float *A, const int LDA, const float *B, const int LDB, const float BETA, float *C, const int LDC) const
 Epetra_BLAS symmetric matrix-matrix multiply function (SSYMM)
 
void SYMM (const char SIDE, const char UPLO, const int M, const int N, const double ALPHA, const double *A, const int LDA, const double *B, const int LDB, const double BETA, double *C, const int LDC) const
 Epetra_BLAS matrix-matrix multiply function (DSYMM)
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const float ALPHA, const float *A, const int LDA, float *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (STRMM)
 
void TRMM (const char SIDE, const char UPLO, const char TRANSA, const char DIAG, const int M, const int N, const double ALPHA, const double *A, const int LDA, double *B, const int LDB) const
 Epetra_BLAS triangular matrix-matrix multiply function (DTRMM)
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const float ALPHA, const float *A, const int LDA, const float BETA, float *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (ssyrk)
 
void SYRK (const char UPLO, const char TRANS, const int N, const int K, const double ALPHA, const double *A, const int LDA, const double BETA, double *C, const int LDC) const
 Eperta_BLAS symetric rank k funtion (dsyrk)
 
- Public Member Functions inherited from Epetra_RowMatrix
virtual ~Epetra_RowMatrix ()
 Destructor.
 
- Public Member Functions inherited from Epetra_Operator
virtual ~Epetra_Operator ()
 Destructor.
 

Additional Inherited Members

- Static Public Member Functions inherited from Epetra_Object
static void SetTracebackMode (int TracebackModeValue)
 Set the value of the Epetra_Object error traceback report mode. More...
 
static int GetTracebackMode ()
 Get the value of the Epetra_Object error report mode.
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting.
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 
- Protected Member Functions inherited from Epetra_CrsMatrix
bool Allocated () const
 
int SetAllocated (bool Flag)
 
double ** Values () const
 
double * All_Values () const
 
double * Values (int LocalRow) const
 
void InitializeDefaults ()
 
int Allocate ()
 
int InsertValues (int LocalRow, int NumEntries, double *Values, int *Indices)
 
int InsertValues (int LocalRow, int NumEntries, const double *Values, const int *Indices)
 
int InsertValues (int LocalRow, int NumEntries, double *Values, long long *Indices)
 
int InsertValues (int LocalRow, int NumEntries, const double *Values, const long long *Indices)
 
int InsertOffsetValues (long long GlobalRow, int NumEntries, double *Values, int *Indices)
 
int InsertOffsetValues (long long GlobalRow, int NumEntries, const double *Values, const int *Indices)
 
int ReplaceOffsetValues (long long GlobalRow, int NumEntries, const double *Values, const int *Indices)
 
int SumIntoOffsetValues (long long GlobalRow, int NumEntries, const double *Values, const int *Indices)
 
void UpdateImportVector (int NumVectors) const
 
void UpdateExportVector (int NumVectors) const
 
void GeneralMV (double *x, double *y) const
 
void GeneralMTV (double *x, double *y) const
 
void GeneralMM (double **X, int LDX, double **Y, int LDY, int NumVectors) const
 
void GeneralMTM (double **X, int LDX, double **Y, int LDY, int NumVectors) const
 
void GeneralSV (bool Upper, bool Trans, bool UnitDiagonal, double *x, double *y) const
 
void GeneralSM (bool Upper, bool Trans, bool UnitDiagonal, double **X, int LDX, double **Y, int LDY, int NumVectors) const
 
void SetStaticGraph (bool Flag)
 
int CheckSizes (const Epetra_SrcDistObject &A)
 Allows the source and target (this) objects to be compared for compatibility, return nonzero if not.
 
int CopyAndPermute (const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)
 Perform ID copies and permutations that are on processor.
 
int CopyAndPermuteCrsMatrix (const Epetra_CrsMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode)
 
int CopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode)
 
int PackAndPrepare (const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)
 Perform any packing or preparation required for call to DoTransfer().
 
int UnpackAndCombine (const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor)
 Perform any unpacking and combining after call to DoTransfer().
 
int SortEntries ()
 Sort column entries, row-by-row, in ascending order.
 
bool Sorted () const
 If SortEntries() has been called, this query returns true, otherwise it returns false.
 
int MergeRedundantEntries ()
 Add entries that have the same column index. Remove redundant entries from list.
 
bool NoRedundancies () const
 If MergeRedundantEntries() has been called, this query returns true, otherwise it returns false.
 
void DeleteMemory ()
 
- Protected Member Functions inherited from Epetra_DistObject
virtual int DoTransfer (const Epetra_SrcDistObject &A, Epetra_CombineMode CombineMode, int NumSameIDs, int NumPermuteIDs, int NumRemoteIDs, int NumExportIDs, int *PermuteToLIDs, int *PermuteFromLIDs, int *RemoteLIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &LenImports, char *&Imports, Epetra_Distributor &Distor, bool DoReverse, const Epetra_OffsetIndex *Indexor)
 Perform actual transfer (redistribution) of data across memory images, using Epetra_Distributor object.
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 
- Protected Attributes inherited from Epetra_CrsMatrix
Epetra_CrsGraph Graph_
 
bool Allocated_
 
bool StaticGraph_
 
bool UseTranspose_
 
bool constructedWithFilledGraph_
 
bool matrixFillCompleteCalled_
 
bool StorageOptimized_
 
double ** Values_
 
int * Values_alloc_lengths_
 
double * All_Values_
 
double NormInf_
 
double NormOne_
 
double NormFrob_
 
int NumMyRows_
 
Epetra_MultiVectorImportVector_
 
Epetra_MultiVectorExportVector_
 
Epetra_DataAccess CV_
 
bool squareFillCompleteCalled_
 
- Protected Attributes inherited from Epetra_DistObject
Epetra_BlockMap Map_
 
const Epetra_CommComm_
 
char * Exports_
 
char * Imports_
 
int LenExports_
 
int LenImports_
 
int * Sizes_
 
- Protected Attributes inherited from Epetra_CompObject
Epetra_FlopsFlopCounter_
 

Detailed Description

Epetra_OskiMatrix: A class for constructing and using OSKI Matrices within Epetra. For information on known issues with OSKI see the detailed description.

OSKI is a high-performance sparse matrix kernel package written by the UC Berkeley Benchmarking and Optimization Group. The Epetra_OskiMatrix class is a lightweight interface to allow Epetra users access to OSKI functionality. This interface includes lightweight conversion of Epetra matrices to OSKI matrices, runtime tuning based on parameter list hints, and access to high-performance computational kernels to perform matrix-vector/matrix multi-vector calculations.

The Epetra_OskiMatrix class provides access to the entire OSKI interface. However, the following features are not fully implemented in OSKI version oski-1.0.1h, and therefore are unsupported in the Epetra_OskiMatrix class:

Constructor & Destructor Documentation

Epetra_OskiMatrix::Epetra_OskiMatrix ( const Epetra_CrsMatrix Source,
const Teuchos::ParameterList &  List 
)

Constructor creates an Epetra_OskiMatrix from an Epetra_CrsMatrix.

Parameters
Source(In) An Epetra_CrsMatrix that is to be wrapped as an Epetra_OskiMatrix.
List(In) Any options or data wanted or needed for the conversion.
Returns
Pointer to an Epetra_OskiMatrix.

Options that can be passed to the List are presented below. They are: "<type> <option name> <default value>: <description of purpose>"

  • bool autotune false: If true, Epetra tries to set as many hints as possible based on its knowledge of the matrix.

string matrixtype general: Other types that can be taken are: uppertri, lowertri, uppersymm, lowersymm, fullsymm, upperherm, lowerherm and fullherm.

  • bool diagstored false: If true, the diagonal entries are not stored in the matrix and are all assumed to be 1.

bool zerobased false: If true, the array is zero based, as in C. Otherwise, it is 1 based, as in Fortran.

  • bool sorted false: If true, all elements in the passed in array are sorted.
  • bool unique false: If true, a value in a column only appears once in each row.
  • bool deepcopy false: If true, when the OSKI matrix is created it will be a deepcopy of the data in the function.

Member Function Documentation

int Epetra_OskiMatrix::ApplyMatrixTransforms ( const char *  Transforms)

Replaces the current data structure of the matrix with the one specified in Transforms.

If a previously tuned copy of the Matrix existed it is now replaced by one specified in Transforms.

Parameters
Transforms(In) A string that holds the transformations to be applied to the matrix. If Transforms is NULL or the empty string then no changes to the data structure are made.
Returns
If the transformation was successful 0 is returned. Otherwise an error code is returned.
int Epetra_OskiMatrix::ExtractDiagonalCopy ( Epetra_Vector Diagonal) const
virtual

Returns a copy of the main diagonal in a user-provided vector.

Parameters
Diagonal- (Out) Extracted main diagonal.
Returns
Integer error code, set to 0 if successful and non-zero if not.
Precondition
Filled()==true
Postcondition
Unchanged.

Reimplemented from Epetra_CrsMatrix.

char* Epetra_OskiMatrix::GetMatrixTransforms ( ) const

Returns a string holding the transformations performed on the matrix when it was tuned.

Returns
Upon success returns a newly-allocated string that stores the transformations applied to the matrix during tuning. NULL is returned upon an error. It is the users responsibility to deallocate the returned string.
int Epetra_OskiMatrix::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 implemented as described in the detailed description.

The vectors x and y can be either Epetra_Vectors or Epetra_OskiVectors. The vector T can be either an Epetra_MultiVector or Epetra_OskiMultiVector. This composed routine is used in power and S-step methods. This routine is not implemented due a bug in the oski-1.01h kernel that makes testing of correctness impossible.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the product.
x(In) The vector the matrix is multiplied by.
y(In/Out) The vector where the calculation result is stored.
T(Out) The multi-vector where the result of each subsequent multiplication this*x ... this^(Power-1)*x is stored.
Power(In) The power to raise the matrix to in the calculation.
Alpha(In) A scalar constant used to scale x.
Beta(In) A scalar constant used to scale y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::MatPowMultiply ( bool  TransA,
const Epetra_Vector x,
Epetra_Vector y,
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 implemented as described in the detailed description.

The vectors x and y can be either Epetra_Vectors or Epetra_OskiVectors. This composed routine is used in power and S-step methods. This routine is not implemented due a bug in the oski-1.01h kernel that makes testing of correctness impossible.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the product.
x(In) The vector the matrix is multiplied by.
y(In/Out) The vector where the calculation result is stored.
Power(In) The power to raise the matrix to in the calculation.
Alpha(In) A scalar constant used to scale x.
Beta(In) A scalar constant used to scale y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::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^TransA*x + Beta*y.

The vectors x, y and t can be either Epetra_Vectors or Epetra_OskiVectors. This composed routine is most commonly used in linear least squares and bidiagonalization methods. The parallel version of y = Alpha*this*this^TransA*x + Beta*y uses calls to the Multiply routine under the hood, as it is not possible to perform both multiplies automatically.

Parameters
ATA(In) If TransA = TRUE then compute this^T*this*x otherwise compute this*this^T*x.
x(In) The vector the matrix is multiplied by.
y(In/Out) The vector where the calculation result is stored.
t(Out) The vector where the result of the this*x is stored if TransA = true and this^T*x is stored otherwise.
Alpha(In) A scalar constant used to scale x.
Beta(In) A scalar constant used to scale y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::MatTransMatMultiply ( bool  ATA,
const Epetra_MultiVector X,
Epetra_MultiVector Y,
Epetra_MultiVector T,
double  Alpha = 1.0,
double  Beta = 0.0 
) const

Performs two matrix multi-vector multiplies of Y = Alpha*this^TransA*this*X + Beta*Y or Y = Alpha*this*this^TransA*X + Beta*Y.

The multi-vectors X, Y and T can be either Epetra_MultiVectors or Epetra_OskiMultiVectors. This composed routine is most commonly used in linear least squares and bidiagonalization methods. The parallel version of Y = Alpha*this*this^TransA*X + Beta*Y uses calls to the Multiply routine under the hood, as it is not possible to perform both multiplies automatically.

Parameters
ATA(In) If TransA = TRUE then compute this^T*this*X otherwise compute this*this^T*X.
X(In) The vector the matrix is multiplied by.
Y(In/Out) The vector where the calculation result is stored.
T(Out) The multi-vector where the result of the this*X is stored if TransA = true and this^T*X is stored otherwise.
Alpha(In) A scalar constant used to scale X.
Beta(In) A scalar constant used to scale Y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::Multiply ( bool  TransA,
const Epetra_Vector x,
Epetra_Vector y 
) const

Performs a matrix vector multiply of y = this^TransA*x.

The vectors x and y can be either Epetra_Vectors or Epetra_OskiVectors.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the product.
x(In) The vector the matrix is multiplied by.
y(Out) The vector where the calculation result is stored.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::Multiply ( bool  TransA,
const Epetra_Vector x,
Epetra_Vector y,
double  Alpha,
double  Beta = 0.0 
) const

Performs a matrix vector multiply of y = Alpha*this^TransA*x + Beta*y.

The vectors x and y can be either Epetra_Vectors or Epetra_OskiVectors.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the product.
x(In) The vector the matrix is multiplied by.
y(In/Out) The vector where the calculation result is stored.
Alpha(In) A scalar constant used to scale x.
Beta(In) A scalar constant used to scale y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Performs a matrix multi-vector multiply of Y = this^TransA*X.

The multi-vectors X and Y can be either Epetra_MultiVectors or Epetra_OskiMultiVectors.

Parameters
TransA(In) If Trans = TRUE then use the transpose of the matrix in computing the product.
X(In) The multi-vector the matrix is multiplied by.
Y(Out) The multi-vector where the calculation result is stored.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged

Reimplemented from Epetra_CrsMatrix.

int Epetra_OskiMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y,
double  Alpha,
double  Beta = 0.0 
) const

Performs a matrix multi-vector multiply of Y = Alpha*this^TransA*X + Beta*Y.

The multi-vectors X and Y can be either Epetra_MultiVectors or Epetra_OskiMultiVectors.

Parameters
TransA(In) If Trans = TRUE then use the transpose of the matrix in computing the product.
X(In) The multi-vector the matrix is multiplied by.
Y(In/Out) The multi-vector where the calculation result is stored.
Alpha(In) A scalar constant used to scale X.
Beta(In) A scalar constant used to scale Y.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::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 + Zeta*z.

The vectors x, y, w and z can be either Epetra_Vectors or Epetra_OskiVectors. This composed routine is most commonly used in bi-conjugate gradient calculations.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the second product.
x(In) A vector the matrix is multiplied by.
y(In/Out) A vector where the calculation result of the first multiply is stored.
w(In) A vector the matrix is multiplied by.
z(In/Out) A vector where the calculation result of the second multiply is stored.
Alpha(In) A scalar constant used to scale x.
Beta(In) A scalar constant used to scale y.
Omega(In) A scalar constant used to scale w.
Zeta(In) A scalar constant used to scale z.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::MultiplyAndMatTransMultiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y,
const Epetra_MultiVector W,
Epetra_MultiVector Z,
double  Alpha = 1.0,
double  Beta = 0.0,
double  Omega = 1.0,
double  Zeta = 0.0 
) const

Performs the two matrix multi-vector multiplies of Y = Alpha*this*X + Beta*Y and Z = Omega*this^TransA*W + Zeta*Z.

The multi-vectors X, Y, W and Z can be either Epetra_MultiVectors or Epetra_OskiMultiVectors. This composed routine is most commonly used in bi-conjugate gradient calculations.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in computing the second product.
X(In) A multi-vector the matrix is multiplied by.
Y(In/Out) A multi-vector where the calculation result of the first multiply is stored.
W(In) A multi-vector the matrix is multiplied by.
Z(In/Out) A multi-vector where the calculation result of the second multiply is stored.
Alpha(In) A scalar constant used to scale X.
Beta(In) A scalar constant used to scale Y.
Omega(In) A scalar constant used to scale W.
Zeta(In) A scalar constant used to scale Z.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::ReplaceDiagonalValues ( const Epetra_OskiVector Diagonal)

Replaces diagonal values of the matrix with those in the user-provided vector.

This routine is meant to allow replacement of { existing} diagonal values. If a diagonal value does not exist for a given row, the corresponding value in the input Epetra_OskiVector will be ignored, and the return code will be set to 1.

Parameters
Diagonal- (In) New values to be placed in the main diagonal.
Returns
Integer error code, set to 0 if successful, set to 1 on the calling processor if one or more diagonal entries not present in matrix. Other error codes can be returned as well, indicating improperly constructed matrices or vectors.
Precondition
Filled()==true
Postcondition
Diagonal values have been replaced with the values of Diagonal.
int Epetra_OskiMatrix::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.

The reason this function could be expensive is its underlying implementation. Both the OSKI and Epetra versions of the matrix must be changed when the matrix has been permuted. When this is the case, a call must be made to the Epetra ReplaceMyValues, and NumEntries calls must be made to a function that changes the OSKI matrix's values one at a time.

Parameters
MyRow(In) Row number (in local coordinates) to put elements.
NumEntries(In) Number of entries.
Values(In) Values to enter.
Indices(In) Local column indices corresponding to the values.
Returns
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, Oski will fail. A positive warning code may be returned, but this should be treated as a fatal error; part of the data will be changed, and OSKI cannot support adding in new data values.
Precondition
IndicesAreLocal()==true
Postcondition
The given Values at the given Indices have been summed into the entries of MyRow.
int Epetra_OskiMatrix::SetHint ( const Teuchos::ParameterList &  List)

Stores the hints in List in the matrix structure.

    \param List (In) A list of hints and options to register along with the matrix
       used for tuning purposes.  The full list is given below.  It may
       be moved to the user guide in the future.
    \return On successful storage of the hint 0 is returned.  On failure an error code
        is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix. For options where multiple arguments can be passed in at once the interface only supports up to 5. This means only 5 block sizes or 5 diaganols can be passed in at once. If you have more changing the code to support your needs should be simple, but every case you use must be enumerated. Of course you can just say there are diagonals and blocks and not pass in specific sizes as wells.

    - bool noblocks: If true, the matrix has no block structure
    - bool singleblocksize: If true, the matrix structure is dominated by blocks of the size of the next two parameters.
          - int row: The number of rows in each block.
          - int col: The number of columns in each block.
    - bool multipleblocksize: If true, the matrix consists of multiple block sizes.  The next 3 parameters describe these and are optional.
      - int blocks: The number of block sizes in the matrix.
      - int row<x>: Where x is the block number, and x goes from 1 to blocks.  This is the number of rows in block x.
      - int col<x>: Where x is the block number, and x goes from 1 to blocks.  This is the number of cols in block x.

    - bool alignedblocks: If true, all blocks are aligned to a grid.
    - bool unalignedblocks: If true, blocks are not aligned to a grid.

    - bool symmetricpattern: If true, the matrix is either symmetric or nearly symmetric.
    - bool nonsymmetricpattern: If true, the matrix has a very unsymmetric pattern.

    - bool randompattern: If true, the matrix's non-zeros are distributed in a random pattern.
    - bool correlatedpattern: If true, the row and column indices for
      non-zeros are highly correlated.

    - bool nodiags : If true, the matrix has little or no diagonal structure.
    - bool diags: If true, the matrix consists of diagonal structure described the next two optional parameters.
          - int numdiags: The number of diagonal sizes known to be present others not listed could be present.
      - int diag<x>: Where x is the diagonal number, and x goes from 1 to numdiags.  This is the size of the diagonal.
int Epetra_OskiMatrix::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 optimize the data structure storage, and the routine to compute the calculation.

    In parallel the routine uses symbolic vectors.  This is done for two reasons.  Doing this
        saves on data allocation and potentially communication overhead.  For a matrix-vector
        routine there should be no advantage to having the actual vector, as its size must be the same
        as a matrix dimension.  For a matrix-multivector routine there could be gains from knowing the
        number of vectors in the multi-vector. However, OSKI does not perform multi-vector optimizations,
        so there is no need to add the overhead.
        \param ATA (In) If ATA = true then this^T*this*x will be computed otherwise this*this^T*x will be.
    \param Alpha (In) A scalar constant used to scale InVec.
    \param InVec (In) The vector the matrix is multiplied by or whether it is a single vector or multi-vector.
    \param Beta  (In) A scalar constant used to scale OutVec.
    \param OutVec (In) The vector where the calculation result is stored or whether it is a single vector or multi-vector.
    \param Intermediate (In) The vector where result of the first product can be stored
       or whether it is a single vector or multi-vector.  If this quantity is NULL
       then the intermediate product is not stored.
    \param NumCalls (In) The number of times the operation is called or the tuning level wanted.
    \param List (In) Used for denoting the use of symbolic vectors for InVec,
       OutVec and Intermediate, along with the level of aggressive tuning if either NumCalls not
       known or to be overridden.  Options are shown below. It should be noted that
       by using these options the associated vector or NumCalls becomes invalid.
    \return Stores the workload hint in the matrix if the operation is valid.  If the
            operation is not valid an error code is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix.

    These replace InVec.
    - bool symminvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symminmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace OutVec.
    - bool symmoutvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmoutmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace Intermediate.
    - bool symmintervec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmintermultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    - bool tune: If true, have OSKI tune moderately rather than using the number of calls passed in.
    - bool tuneaggressive: If true, have OSKI tune aggressively rather than using the number of calls passed in.
int Epetra_OskiMatrix::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 structure storage, and the routine to compute the calculation.

    In parallel the routine uses symbolic vectors.  This is done for two reasons.  Doing this
        saves on data allocation and potentially communication overhead.  For a matrix-vector
        routine there should be no advantage to having the actual vector, as its size must be the same
        as a matrix dimension.  For a matrix-multivector routine there could be gains from knowing the
        number of vectors in the multi-vector. However, OSKI does not perform multi-vector optimizations,
        so there is no need to add the overhead.
        \param Trans (In) If Trans = true then the transpose of the matrix will be used in
           computing the product.
    \param Alpha (In) A scalar constant used to scale InVec.
    \param InVec (In) The vector the matrix is multiplied by or whether it is a single vector or multi-vector.
    \param Beta  (In) A scalar constant used to scale OutVec.
    \param OutVec (In) The vector where the calculation result is stored or whether it is a single vector or multi-vector.
    \param NumCalls (In) The number of times the operation is called or the tuning level wanted.
    \param List (In) Used for denoting the use of symbolic vectors for both InVec
       and OutVec, as well as for level of aggressive tuning if either NumCalls not
       known or to be overridden.  Options are shown below. It should be noted that
       by using these options the associated vector or NumCalls becomes invalid.
    \return Stores the workload hint in the matrix if the operation is valid.  If the
            operation is not valid an error code is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix.

    These replace InVec.
    - bool symminvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symminmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace OutVec.
    - bool symmoutvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmoutmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    - bool tune: If true, have OSKI tune moderately rather than using the number of calls passed in.
    - bool tuneaggressive: If true, have OSKI tune aggressively rather than using the number of calls passed in.
int Epetra_OskiMatrix::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 structure storage, and the routine to compute the calculation.

    In parallel the routine uses symbolic vectors.  This is done for two reasons.  Doing this
        saves on data allocation and potentially communication overhead.  For a matrix-vector
        routine there should be no advantage to having the actual vector, as its size must be the same
        as a matrix dimension.  For a matrix-multivector routine there could be gains from knowing the
        number of vectors in the multi-vector. However, OSKI does not perform multi-vector optimizations,
        so there is no need to add the overhead.

        \param Trans (In) If Trans = true then the transpose of the matrix will be used in
           computing the product.
    \param Alpha (In) A scalar constant used to scale InVec.
    \param InVec (In) The vector the matrix is multiplied by or whether it is a single vector or multi-vector.
    \param Beta  (In) A scalar constant used to scale OutVec.
    \param OutVec (In) The vector where the calculation result is stored or whether it is a single vector or multi-vector.
    \param Omega (In) A scalar constant used to scale InVec2.
    \param InVec2 (In) The vector the matrix is multiplied by or whether it is a single vector or multi-vector.
    \param Zeta  (In) A scalar constant used to scale OutVec2.
    \param OutVec2 (In) The vector where the calculation result is stored or whether it is a single vector or multi-vector.
    \param NumCalls (In) The number of times the operation is called or the tuning level wanted.
    \param List (In) Used for denoting the use of symbolic vectors for both InVec
       and OutVec, as well as for level of aggressive tuning if either NumCalls not
       known or to be overridden.  Options are shown below. It should be noted that
       by using these options the associated vector or NumCalls becomes invalid.
    \return Stores the workload hint in the matrix if the operation is valid.  If the
            operation is not valid an error code is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix.

    These replace InVec.
    - bool symminvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symminmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace OutVec.
    - bool symmoutvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmoutmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    These replace InVec2.
    - bool symminvec2: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symminmultivec2: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace OutVec2.
    - bool symmoutvec2: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmoutmultivec2: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    - bool tune: If true, have OSKI tune moderately rather than using the number of calls passed in.
    - bool tuneaggressive: If true, have OSKI tune aggressively rather than using the number of calls passed in.
int Epetra_OskiMatrix::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 optimize the data structure storage and the routine to compute the calculation.

    In parallel the routine uses symbolic vectors.  This is done for two reasons.  Doing this
        saves on data allocation and potentially communication overhead.  For a matrix-vector
        routine there should be no advantage to having the actual vector, as its size must be the same
        as a matrix dimension.  For a matrix-multivector routine there could be gains from knowing the
        number of vectors in the multi-vector. However, OSKI does not perform multi-vector optimizations,
        so there is no need to add the overhead.
        \param Trans (In) If Trans = true then the transpose of the matrix will be used in
           computing the product.
    \param Alpha (In) A scalar constant used to scale InVec.
    \param InVec (In) The vector the matrix is multiplied by or whether it is a single vector or multi-vector.
    \param Beta  (In) A scalar constant used to scale OutVec.
    \param OutVec (In) The vector where the calculation result is stored or whether it is a single vector or multi-vector.
    \param Intermediate (In) The multi-vector where result of the first product can be stored
       or whether it is a single vector or multi-vector.  If this quantity is NULL
       then the intermediate product is not stored.
    \param Power (In) The power to raise the matrix to in the calculation.
    \param NumCalls (In) The number of times the operation is called or the tuning level wanted.
    \param List (In) Used for denoting the use of symbolic vectors for both InVec
       and OutVec, as well as for level of aggressive tuning if either NumCalls not
       known or to be overridden.  Options are shown below. It should be noted that
       by using these options the associated vector or NumCalls becomes invalid.
    \return Stores the workload hint in the matrix if the operation is valid.  If the
            operation is not valid an error code is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix.

    These replace InVec.
    - bool symminvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symminmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        These replace OutVec.
    - bool symmoutvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmoutmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

        This replaces Intermediate.
    - bool symmintermultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    - bool tune: If true, have OSKI tune moderately rather than using the number of calls passed in.
    - bool tuneaggressive: If true, have OSKI tune aggressively rather than using the number of calls passed in.
int Epetra_OskiMatrix::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 storage, and the routine to compute the calculation.

    In parallel the routine uses symbolic vectors.  This is done for two reasons.  Doing this
        saves on data allocation and potentially communication overhead.  For a matrix-vector
        routine there should be no advantage to having the actual vector, as its size must be the same
        as a matrix dimension.  For a matrix-multivector routine there could be gains from knowing the
        number of vectors in the multi-vector. However, OSKI does not perform multi-vector optimizations,
        so there is no need to add the overhead.

        \param Trans (In) If Trans = true then the transpose of the matrix will be used in
           computing the product.
    \param Alpha (In) A scalar constant used to scale InVec.
    \param Vector (In) The vector being used in the solve and to store the solution.
    \param NumCalls (In) The number of times the operation is called or the tuning level wanted.
    \param List (In) Used for denoting the use of a symbolic vectors, as well as for
       level of aggressive tuning if either NumCalls not
       known or to be overridden.  Options are shown below. It should be noted that
       by using these options the associated vector or NumCalls becomes invalid.
    \return Stores the workload hint in the matrix if the operation is valid.  If the
            operation is not valid an error code is returned.

        Options that can be passed to the List are presented below.  They are: "<type> <option name>: <description of purpose>".  The

available hints are grouped by section, and only one hint from each section can be true for a given matrix.

    These replace Vector.
    - bool symmvec: If true, use a symbolic vector rather than the vector passed in for tuning purposes.
    - bool symmmultivec: If true, use a symbolic multi-vector rather than the multi-vector passed in for tuning purposes.

    - bool tune: If true, have OSKI tune moderately rather than using the number of calls passed in.
    - bool tuneaggressive: If true, have OSKI tune aggressively rather than using the number of calls passed in.
int Epetra_OskiMatrix::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.

The vector x can be either be an Epetra_Vector or Epetra_OskiVector. The OskiMatrix must already be triangular, and the UnitDiagonal setting associated with it will be used.

Parameters
Upper(In) This parameter is ignored, and is here only to match the Epetra_CrsMatrix::Solve syntax.
TransA(In) If TransA = TRUE then use the transpose of the matrix in solving the equations.
UnitDiagonal(In) This parameter is ignored only here to match the Epetra_CrsMatrix::Solve syntax.
x(In) The vector solved against.
y(Out) The solution vector.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::Solve ( bool  TransA,
const Epetra_Vector x,
Epetra_Vector y,
double  Alpha = 1.0 
) const

Performs a triangular solve of y = Alpha*(this^TransA)^-1*x where this is a triangular matrix.

The vector x can be either be an Epetra_Vector or Epetra_OskiVector.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in solving the equations.
x(In) The vector solved against.
y(Out) The solution vector.
Alpha(In) A scalar constant used to scale x.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::Solve ( bool  Upper,
bool  TransA,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Performs a triangular solve of Y = (this^TransA)^-1*X where this is a triangular matrix.

The vector X can be either be an Epetra_MultiVector or Epetra_OskiMultiVector. The OskiMatrix must already be triangular, and the UnitDiagonal setting associated with it will be used.

Parameters
Upper(In) This parameter is ignored only here to match the Epetra_CrsMatrix::Solve syntax.
TransA(In) If TransA = TRUE then use the transpose of the matrix in solving the equations.
UnitDiagonal(In) This parameter is ignored only here to match the Epetra_CrsMatrix::Solve syntax.
X(In) The multi-vector solved against.
Y(Out) The solution multi-vector.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged

Reimplemented from Epetra_CrsMatrix.

int Epetra_OskiMatrix::Solve ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y,
double  Alpha = 1.0 
) const

Performs a triangular solve of Y = Alpha*(this^TransA)^-1*X where this is a triangular matrix.

The vector X can be either be an Epetra_MultiVector or Epetra_OskiMultiVector.

Parameters
TransA(In) If TransA = TRUE then use the transpose of the matrix in solving the equations.
X(In) The multi-vector solved against.
Y(Out) The solution multi-vector.
Alpha(In) A scalar constant used to scale X.
Returns
Integer error code, set to 0 if successful.
Precondition
Filled()==true
Postcondition
Unchanged
int Epetra_OskiMatrix::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 be expensive.

The reason this function could be expensive is its underlying implementation. Both the OSKI and Epetra versions of the Matrix must be changed when the matrix has been permuted. When this is the case, a call must be made to the Epetra SumIntoMyValues, and NumEntries calls must be made to a function that changes the OSKI matrix's values one at a time.

Parameters
MyRow- (In) Row number (in local coordinates) to put elements.
NumEntries- (In) Number of entries.
Values- (In) Values to enter.
Indices- (In) Local column indices corresponding to values.
Returns
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code may be returned. This should be treated as a fatal error, as part of the data will be changed, and OSKI cannot support adding in new data values.
Precondition
IndicesAreLocal()==true
Postcondition
The given Values at the given Indices have been summed into the entries of MyRow.
int Epetra_OskiMatrix::TuneMatrix ( )

Tunes the matrix multiply if its deemed profitable.

The routine tunes based upon user provided hints if given. If hints are not given the tuning is performed based on expected future workload for the calculation.

Returns
On success returns a non-negative status code of the transformations performed. On failure an error code is returned.

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