#include <Epetra_CrsMatrix.h>
Inheritance diagram for Epetra_CrsMatrix:


Public Member Functions | |
Constructors/Destructor | |
| 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. | |
| Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int NumEntriesPerRow, bool StaticProfile=false) | |
| Epetra_CrsMatrix constructor with fixed number of indices per row. | |
| 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. | |
| 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. | |
| Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_CrsGraph &Graph) | |
| Construct a matrix using an existing Epetra_CrsGraph object. | |
| Epetra_CrsMatrix (const Epetra_CrsMatrix &Matrix) | |
| Copy constructor. | |
| virtual | ~Epetra_CrsMatrix () |
| Epetra_CrsMatrix Destructor. | |
Insertion/Replace/SumInto methods | |
| Epetra_CrsMatrix & | operator= (const Epetra_CrsMatrix &src) |
| Assignment operator. | |
| int | PutScalar (double ScalarConstant) |
| Initialize all values in the matrix with constant value. | |
| int | Scale (double ScalarConstant) |
| Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A). | |
| virtual int | InsertGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| Insert a list of elements in a given global row of the matrix. | |
| virtual int | ReplaceGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| Replace specified existing values with this list of entries for a given global row of the matrix. | |
| virtual int | SumIntoGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| Add this list of entries to existing values for a given global row of the matrix. | |
| int | InsertMyValues (int MyRow, int NumEntries, double *Values, int *Indices) |
| Insert a list of elements in a given local row of the matrix. | |
| 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. | |
| 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. | |
| int | ReplaceDiagonalValues (const Epetra_Vector &Diagonal) |
| Replaces diagonal values of the matrix with those in the user-provided vector. | |
Transformation methods | |
| 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. | |
| int | MakeDataContiguous () |
| Eliminates memory that is used for construction. Make consecutive row index sections contiguous. | |
Extraction methods | |
| 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. | |
| 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. | |
| int | ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values) const |
| Returns a copy of the specified global row values in user-provided array. | |
| int | ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values) const |
| Returns a copy of the specified local row values in user-provided array. | |
| int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const |
| Returns a copy of the main diagonal in a user-provided vector. | |
| 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. | |
| 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. | |
| int | ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values) const |
| Returns a view of the specified global row values via pointers to internal data. | |
| int | ExtractMyRowView (int MyRow, int &NumEntries, double *&Values) const |
| Returns a view of the specified local row values via pointers to internal data. | |
Computational methods | |
| 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. | |
| int | Multiply1 (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const |
| int | Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y. | |
| 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. | |
| int | Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
| Returns the result of a local solve using the Epetra_CrsMatrix a Epetra_MultiVector X in Y. | |
| 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. | |
| 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. | |
| int | LeftScale (const Epetra_Vector &x) |
| Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x. | |
| 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. | |
| int | InvColMaxs (Epetra_Vector &x) const |
| Computes the max of absolute values of the columns of the Epetra_CrsMatrix, results returned in x. | |
| int | RightScale (const Epetra_Vector &x) |
| Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x. | |
Matrix Properties Query Methods | |
| 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. | |
Atribute access functions | |
| 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. | |
| int | NumGlobalRows () const |
| Returns the number of global matrix rows. | |
| int | NumGlobalCols () const |
| Returns the number of global matrix columns. | |
| int | NumGlobalDiagonals () const |
| Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. | |
| 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. | |
| int | NumMyDiagonals () const |
| Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. | |
| int | NumGlobalEntries (int 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. | |
| int | GlobalMaxNumEntries () const |
| Returns the maximum number of nonzero entries across all rows on all processors. | |
| 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. | |
| bool | StaticGraph () |
| Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable. | |
| const Epetra_CrsGraph & | Graph () const |
| Returns a reference to the Epetra_CrsGraph object associated with this matrix. | |
| const Epetra_Map & | RowMap () 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. | |
| bool | HaveColMap () const |
| Returns true if we have a well-defined ColMap, and returns false otherwise. | |
| int | ReplaceColMap (const Epetra_BlockMap &newmap) |
| Replaces the current ColMap with the user-specified map object. | |
| const Epetra_Map & | ColMap () const |
| Returns the Epetra_Map object that describes the set of column-indices that appear in each processor's locally owned matrix rows. | |
| const Epetra_Map & | DomainMap () const |
| Returns the Epetra_Map object associated with the domain of this matrix operator. | |
| const Epetra_Map & | RangeMap () const |
| Returns the Epetra_Map object associated with the range of this matrix operator. | |
| const Epetra_Import * | Importer () const |
| Returns the Epetra_Import object that contains the import operations for distributed operations. | |
| const Epetra_Export * | Exporter () const |
| Returns the Epetra_Export object that contains the export operations for distributed operations. | |
| const Epetra_Comm & | Comm () const |
| Returns a pointer to the Epetra_Comm communicator associated with this matrix. | |
Local/Global ID methods | |
| 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 | 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. | |
| 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. | |
| 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. | |
| 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 | 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. | |
| bool | MyLCID (int LCID_in) const |
| Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false. | |
| bool | MyGlobalRow (int GID) const |
| Returns true of GID is owned by the calling processor, otherwise it returns false. | |
I/O Methods | |
| virtual void | Print (ostream &os) const |
| Print method. | |
Additional methods required to support the Epetra_Operator interface | |
| 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. | |
| 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. | |
| 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. | |
| bool | HasNormInf () const |
| Returns true because this class can compute an Inf-norm. | |
| bool | UseTranspose () const |
| Returns the current UseTranspose setting. | |
| const Epetra_Map & | OperatorDomainMap () const |
| Returns the Epetra_Map object associated with the domain of this matrix operator. | |
| const Epetra_Map & | OperatorRangeMap () const |
| Returns the Epetra_Map object associated with the range of this matrix operator. | |
Additional methods required to implement Epetra_RowMatrix interface | |
| int | NumMyRowEntries (int MyRow, int &NumEntries) const |
| Return the current number of values stored for the specified local row. | |
| const Epetra_Map & | RowMatrixRowMap () const |
| Returns the Epetra_Map object associated with the rows of this matrix. | |
| const Epetra_Map & | RowMatrixColMap () const |
| Returns the Epetra_Map object associated with columns of this matrix. | |
| const Epetra_Import * | RowMatrixImporter () const |
| Returns the Epetra_Import object that contains the import operations for distributed operations. | |
Inlined Operator Methods | |
| double * | operator[] (int Loc) |
| Inlined bracket operator for fast access to data. (Const and Non-const versions). | |
| double * | operator[] (int Loc) const |
Expert-only methods: These methods are intended for experts only and have some risk of changing in the future, since they rely on underlying data structure assumptions | |
| int | ExtractCrsDataPointers (int *&IndexOffset, int *&Indices, double *&Values_in) const |
| Returns internal data pointers associated with Crs matrix format. | |
| int | SortGhostsAssociatedWithEachProcessor (bool Flag) |
| Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order. | |
Deprecated methods: These methods still work, but will be removed in a future version | |
| const Epetra_Map & | ImportMap () 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. | |
Protected Member Functions | |
| 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 | InsertOffsetValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| int | ReplaceOffsetValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| int | SumIntoOffsetValues (int GlobalRow, int NumEntries, double *Values, 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) |
| 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) |
| int | CopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor) |
| 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 Attributes | |
| Epetra_CrsGraph | Graph_ |
| bool | Allocated_ |
| bool | StaticGraph_ |
| bool | UseTranspose_ |
| bool | constructedWithFilledGraph_ |
| bool | matrixFillCompleteCalled_ |
| bool | StorageOptimized_ |
| double ** | Values_ |
| double * | All_Values_ |
| double | NormInf_ |
| double | NormOne_ |
| double | NormFrob_ |
| int | NumMyRows_ |
| Epetra_MultiVector * | ImportVector_ |
| Epetra_MultiVector * | ExportVector_ |
| Epetra_DataAccess | CV_ |
| bool | squareFillCompleteCalled_ |
The Epetra_CrsMatrix class is a sparse compressed row matrix object. This matrix can be used in a parallel setting, with data distribution described by Epetra_Map attributes. The structure or graph of the matrix is defined by an Epetra_CrsGraph attribute.
In addition to coefficient access, the primary operations provided by Epetra_CrsMatrix are matrix times vector and matrix times multi-vector multiplication.
Epetra_CrsMatrix matrices can be square or rectangular.
Creating and filling Epetra_CrsMatrix objects
Constructing Epetra_CrsMatrix objects is a multi-step process. The basic steps are as follows:
Note that, even after a matrix is constructed (FillComplete has been called), it is possible to update existing matrix entries. It is not possible to create new entries.
Epetra_Map attributes
Epetra_CrsMatrix objects have four Epetra_Map attributes, which are held by the Epetra_CrsGraph attribute.
The Epetra_Map attributes can be obtained via these accessor methods:
It is important to note that while the row-map and the range-map are often the same, the column-map and the domain-map are almost never the same. The set of entries in a distributed column-map almost always form overlapping sets, with entries being associated with more than one processor. A domain-map, on the other hand, must be a 1-to-1 map, with entries being associated with only a single processor.
Local versus Global Indices
Epetra_CrsMatrix has query functions IndicesAreLocal() and IndicesAreGlobal(), which are used to determine whether the underlying Epetra_CrsGraph attribute's column-indices have been transformed into a local index space or not. (This transformation occurs when the method Epetra_CrsGraph::FillComplete() is called, which happens when the method Epetra_CrsMatrix::FillComplete() is called.) The state of the indices in the graph determines the behavior of many Epetra_CrsMatrix methods. If an Epetra_CrsMatrix instance is constructed using one of the constructors that does not accept a pre-existing Epetra_CrsGraph object, then an Epetra_CrsGraph attribute is created internally and its indices remain untransformed (IndicesAreGlobal()==true) until Epetra_CrsMatrix::FillComplete() is called. The query function Epetra_CrsMatrix::Filled() returns true if Epetra_CrsMatrix::FillComplete() has been called.
Note the following method characteristics:
Most methods have preconditions documented, check documentation for specific methods not mentioned here.
Counting Floating Point Operations
Each Epetra_CrsMatrix object keeps track of the number of serial floating point operations performed using the specified object as the this argument to the function. The Flops() function returns this number as a double precision number. Using this information, in conjunction with the Epetra_Time class, one can get accurate parallel performance numbers. The ResetFlops() function resets the floating point counter.
| Epetra_CrsMatrix::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.
Creates a Epetra_CrsMatrix object and allocates storage.
| CV | - (In) An Epetra_DataAccess enumerated type set to Copy or View. | |
| RowMap | - (In) An Epetra_Map defining the numbering and distribution of matrix rows. | |
| NumEntriesPerRow | - (In) An integer array of length NumRows such that NumEntriesPerRow[i] indicates the (approximate if StaticProfile=false) number of entries in the ith row. | |
| StaticProfile | - (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
| Epetra_CrsMatrix::Epetra_CrsMatrix | ( | Epetra_DataAccess | CV, | |
| const Epetra_Map & | RowMap, | |||
| int | NumEntriesPerRow, | |||
| bool | StaticProfile = false | |||
| ) |
Epetra_CrsMatrix constructor with fixed number of indices per row.
Creates a Epetra_CrsMatrix object and allocates storage.
| CV | - (In) An Epetra_DataAccess enumerated type set to Copy or View. | |
| RowMap | - (In) An Epetra_Map defining the numbering and distribution of matrix rows. | |
| NumEntriesPerRow | - (In) An integer that indicates the (approximate) number of entries in the each row. Note that it is possible to use 0 for this value and let fill occur during the insertion phase. | |
| StaticProfile | - (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
| Epetra_CrsMatrix::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.
Creates a Epetra_CrsMatrix object and allocates storage.
| CV | - (In) An Epetra_DataAccess enumerated type set to Copy or View. | |
| RowMap | - (In) An Epetra_Map defining the numbering and distribution of matrix rows. | |
| ColMap | - (In) An Epetra_Map defining the set of column-indices that appear in each processor's locally owned matrix rows. | |
| NumEntriesPerRow | - (In) An integer array of length NumRows such that NumEntriesPerRow[i] indicates the (approximate if StaticProfile=false) number of entries in the ith row. | |
| StaticProfile | - (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
| Epetra_CrsMatrix::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.
Creates a Epetra_CrsMatrix object and allocates storage.
| CV | - (In) An Epetra_DataAccess enumerated type set to Copy or View. | |
| RowMap | - (In) An Epetra_Map defining the numbering and distribution of matrix rows. | |
| ColMap | - (In) An Epetra_Map defining the set of column-indices that appear in each processor's locally owned matrix rows. | |
| NumEntriesPerRow | - (In) An integer that indicates the (approximate if StaticProfile=false) number of entries in the each row. Note that it is possible to use 0 for this value and let fill occur during the insertion phase. | |
| StaticProfile | - (In) Optional argument that indicates whether or not NumIndicesPerRow should be interpreted as an exact count of nonzeros, or should be used as an approximation. By default this value is false, allowing the profile to be determined dynamically. If the user sets it to true, then the memory allocation for the Epetra_CrsGraph object will be done in one large block, saving on memory fragmentation and generally improving the performance of matrix multiplication and solve kernels. |
| Epetra_CrsMatrix::Epetra_CrsMatrix | ( | Epetra_DataAccess | CV, | |
| const Epetra_CrsGraph & | Graph | |||
| ) |
Construct a matrix using an existing Epetra_CrsGraph object.
Allows the nonzero structure from another matrix, or a structure that was constructed independently, to be used for this matrix.
| CV | - (In) An Epetra_DataAccess enumerated type set to Copy or View. | |
| Graph | - (In) A Epetra_CrsGraph object, constructed directly or extracted from another Epetra matrix object. |
| int Epetra_CrsMatrix::Apply | ( | const Epetra_MultiVector & | X, | |
| Epetra_MultiVector & | Y | |||
| ) | const [inline, virtual] |
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
| X | - (In) An Epetra_MultiVector of dimension NumVectors to multiply with matrix. | |
| Y | -(Out) An Epetra_MultiVector of dimension NumVectors containing result. |
Implements Epetra_Operator.
| int Epetra_CrsMatrix::ApplyInverse | ( | const Epetra_MultiVector & | X, | |
| Epetra_MultiVector & | Y | |||
| ) | const [inline, virtual] |
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
In this implementation, we use several existing attributes to determine how virtual method ApplyInverse() should call the concrete method Solve(). We pass in the UpperTriangular(), the Epetra_CrsMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.
| X | - (In) An Epetra_MultiVector of dimension NumVectors to solve for. | |
| Y | - (Out) An Epetra_MultiVector of dimension NumVectors containing result. |
Implements Epetra_Operator.
| const Epetra_Map& Epetra_CrsMatrix::ColMap | ( | ) | const [inline] |
Returns the Epetra_Map object that describes the set of column-indices that appear in each processor's locally owned matrix rows.
Note that if the matrix was constructed with only a row-map, then until FillComplete() is called, this method returns a column-map that is a copy of the row-map. That 'initial' column-map is replaced with a computed column-map (that contains the set of column-indices appearing in each processor's local portion of the matrix) when FillComplete() is called.
| const Epetra_Map& Epetra_CrsMatrix::DomainMap | ( | ) | const [inline] |
Returns the Epetra_Map object associated with the domain of this matrix operator.
| int Epetra_CrsMatrix::ExtractCrsDataPointers | ( | int *& | IndexOffset, | |
| int *& | Indices, | |||
| double *& | Values_in | |||
| ) | const [inline] |
Returns internal data pointers associated with Crs matrix format.
Returns data pointers to facilitate optimized code within external packages.
| IndexOffset | - (Out) Extracted array of indices into Values[] and Indices[]. Local row k is stored in Values[IndexOffset[k]:IndexOffset[k+1]-1] and Indices[IndexOffset[k]:IndexOffset[k+1]-1]. | |
| Values | - (Out) Extracted values for all local rows. | |
| Indices | - (Out) Extracted local column indices for the corresponding values. |
| int Epetra_CrsMatrix::ExtractDiagonalCopy | ( | Epetra_Vector & | Diagonal | ) | const [virtual] |
Returns a copy of the main diagonal in a user-provided vector.
| Diagonal | - (Out) Extracted main diagonal. |
Implements Epetra_RowMatrix.
Reimplemented in Epetra_OskiMatrix.
| int Epetra_CrsMatrix::ExtractGlobalRowCopy | ( | int | GlobalRow, | |
| int | Length, | |||
| int & | NumEntries, | |||
| double * | Values | |||
| ) | const |
Returns a copy of the specified global row values in user-provided array.
| GlobalRow | - (In) Global row to extract. | |
| Length | - (In) Length of Values. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. |
| int Epetra_CrsMatrix::ExtractGlobalRowCopy | ( | int | GlobalRow, | |
| int | Length, | |||
| int & | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) | const |
Returns a copy of the specified global row in user-provided arrays.
| GlobalRow | - (In) Global row to extract. | |
| ILength | - (In) Length of Values and Indices. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. | |
| Indices | - (Out) Extracted global column indices for the corresponding values. |
| int Epetra_CrsMatrix::ExtractGlobalRowView | ( | int | GlobalRow, | |
| int & | NumEntries, | |||
| double *& | Values | |||
| ) | const |
Returns a view of the specified global row values via pointers to internal data.
| GlobalRow | - (In) Global row to extract. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. |
| int Epetra_CrsMatrix::ExtractGlobalRowView | ( | int | GlobalRow, | |
| int & | NumEntries, | |||
| double *& | Values, | |||
| int *& | Indices | |||
| ) | const |
Returns a view of the specified global row values via pointers to internal data.
| GlobalRow | - (In) Global row to view. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. | |
| Indices | - (Out) Extracted global column indices for the corresponding values. |
| int Epetra_CrsMatrix::ExtractMyRowCopy | ( | int | MyRow, | |
| int | Length, | |||
| int & | NumEntries, | |||
| double * | Values | |||
| ) | const |
Returns a copy of the specified local row values in user-provided array.
| MyRow | - (In) Local row to extract. | |
| Length | - (In) Length of Values. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. |
| int Epetra_CrsMatrix::ExtractMyRowCopy | ( | int | MyRow, | |
| int | Length, | |||
| int & | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) | const [virtual] |
Returns a copy of the specified local row in user-provided arrays.
| MyRow | - (In) Local row to extract. | |
| Length | - (In) Length of Values and Indices. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. | |
| Indices | - (Out) Extracted local column indices for the corresponding values. |
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::ExtractMyRowView | ( | int | MyRow, | |
| int & | NumEntries, | |||
| double *& | Values | |||
| ) | const |
Returns a view of the specified local row values via pointers to internal data.
| MyRow | - (In) Local row to extract. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. |
| int Epetra_CrsMatrix::ExtractMyRowView | ( | int | MyRow, | |
| int & | NumEntries, | |||
| double *& | Values, | |||
| int *& | Indices | |||
| ) | const |
Returns a view of the specified local row values via pointers to internal data.
| MyRow | - (In) Local row to view. | |
| NumEntries | - (Out) Number of nonzero entries extracted. | |
| Values | - (Out) Extracted values for this row. | |
| Indices | - (Out) Extracted local column indices for the corresponding values. |
| int Epetra_CrsMatrix::GCID | ( | int | LCID_in | ) | const [inline] |
Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column.
| int Epetra_CrsMatrix::GlobalMaxNumEntries | ( | ) | const [inline] |
Returns the maximum number of nonzero entries across all rows on all processors.
| bool Epetra_CrsMatrix::HaveColMap | ( | ) | const [inline] |
Returns true if we have a well-defined ColMap, and returns false otherwise.
| virtual int Epetra_CrsMatrix::InsertGlobalValues | ( | int | GlobalRow, | |
| int | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) | [virtual] |
Insert a list of elements in a given global row of the matrix.
This method is used to construct a matrix for the first time. It cannot be used if the matrix structure has already been fixed (via a call to FillComplete()). If multiple values are inserted for the same matrix entry, the values are initially stored separately, so memory use will grow as a result. However, when FillComplete is called the values will be summed together and the additional memory will be released.
For example, if the values 2.0, 3.0 and 4.0 are all inserted in Row 1, Column 2, extra storage is used to store each of the three values separately. In this way, the insert process does not require any searching and can be faster. However, when FillComplete() is called, the values will be summed together to equal 9.0 and only a single entry will remain in the matrix for Row 1, Column 2.
| GlobalRow | - (In) Row number (in global coordinates) to put elements. | |
| NumEntries | - (In) Number of entries. | |
| Values | - (In) Values to enter. | |
| Indices | - (In) Global column indices corresponding to values. |
Reimplemented in Epetra_FECrsMatrix.
| int Epetra_CrsMatrix::InsertMyValues | ( | int | MyRow, | |
| int | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) |
Insert a list of elements in a given local row of the matrix.
| 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. |
| int Epetra_CrsMatrix::InvColMaxs | ( | Epetra_Vector & | x | ) | const |
Computes the max of absolute values of the columns of the Epetra_CrsMatrix, results returned in x.
The vector x will return such that x[j] will contain the inverse of max of the absolute values of the entries in the jth row of the this matrix.
| x | - (Out) An Epetra_Vector containing the column maxs of the this matrix. |
| int Epetra_CrsMatrix::InvColSums | ( | Epetra_Vector & | x | ) | const [virtual] |
Computes the inverse of the sum of absolute values of the columns of the Epetra_CrsMatrix, results returned in x.
The vector x will return such that x[j] will contain the inverse of the sum of the absolute values of the entries in the jth column of the this matrix. Using the resulting vector from this function as input to RightScale() will make the one norm of the resulting matrix exactly 1.
| x | - (Out) An Epetra_Vector containing the column sums of the this matrix. |
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::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.
The vector x will return such that x[i] will contain the inverse of max of the absolute values of the entries in the ith row of the this matrix.
| x | - (Out) An Epetra_Vector containing the inverse of the row maxs of the this matrix. |
| int Epetra_CrsMatrix::InvRowSums | ( | Epetra_Vector & | x | ) | const [virtual] |
Computes the inverse of the sum of absolute values of the rows of the Epetra_CrsMatrix, results returned in x.
The vector x will return such that x[i] will contain the inverse of the sum of the absolute values of the entries in the ith row of the this matrix. Using the resulting vector from this function as input to LeftScale() will make the infinity norm of the resulting matrix exactly 1.
| x | - (Out) An Epetra_Vector containing the inverse of the row sums of the this matrix. |
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::LCID | ( | int | GCID_in | ) | const [inline] |
Returns the local column index for given global column index, returns -1 if no local column for this global column.
| int Epetra_CrsMatrix::LeftScale | ( | const Epetra_Vector & | x | ) | [virtual] |
Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x.
The this matrix will be scaled such that A(i,j) = x(i)*A(i,j) where i denotes the row number of A and j denotes the column number of A.
| x | - (In) An Epetra_Vector to scale with. |
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::MaxNumEntries | ( | ) | const [inline, virtual] |
Returns the maximum number of nonzero entries across all rows on this processor.
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::Multiply | ( | bool | TransA, | |
| const Epetra_MultiVector & | X, | |||
| Epetra_MultiVector & | Y | |||
| ) | const [virtual] |
Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y.
| TransA | - (In) If true, multiply by the transpose of matrix, otherwise just use matrix. | |
| X | - (In) An Epetra_MultiVector of dimension NumVectors to multiply with matrix. | |
| Y | - (Out) An Epetra_MultiVector of dimension NumVectorscontaining result. |
Implements Epetra_RowMatrix.
Reimplemented in Epetra_OskiMatrix.
| int Epetra_CrsMatrix::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.
| TransA | - (In) If true, multiply by the transpose of matrix, otherwise just use matrix. | |
| x | - (In) An Epetra_Vector to multiply by. | |
| y | - (Out) An Epetra_Vector containing result. |
Reimplemented in Epetra_OskiMatrix.
| bool Epetra_CrsMatrix::MyGCID | ( | int | GCID_in | ) | const [inline] |
Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false.
| bool Epetra_CrsMatrix::MyLCID | ( | int | LCID_in | ) | const [inline] |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
| int Epetra_CrsMatrix::NumMyCols | ( | ) | const [inline, virtual] |
Returns the number of entries in the set of column-indices that appear on this processor.
The set of column-indices that appear on this processor is the union of column-indices that appear in all local rows. The size of this set isn't available until FillComplete() has been called.
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::NumMyDiagonals | ( | ) | const [inline, virtual] |
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::NumMyRowEntries | ( | int | MyRow, | |
| int & | NumEntries | |||
| ) | const [virtual] |
Return the current number of values stored for the specified local row.
Similar to NumMyEntries() except NumEntries is returned as an argument and error checking is done on the input value MyRow.
| MyRow | - (In) Local row. | |
| NumEntries | - (Out) Number of nonzero values. |
Implements Epetra_RowMatrix.
| double* Epetra_CrsMatrix::operator[] | ( | int | Loc | ) | [inline] |
Inlined bracket operator for fast access to data. (Const and Non-const versions).
No error checking and dangerous for optimization purposes.
| Loc | - (In) Local row. |
| int Epetra_CrsMatrix::OptimizeStorage | ( | ) |
Make consecutive row index sections contiguous, minimize internal storage used for constructing graph.
After construction and during initialization (when values are being added), the matrix coefficients for each row are managed as separate segments of memory. This method moves the values for all rows into one large contiguous array and eliminates internal storage that is not needed after matrix construction. Calling this method can have a significant impact on memory costs and machine performance.
If this object was constructed in View mode then this method can't make non-contiguous values contiguous and will return a warning code of 1 if the viewed data isn't already contiguous.
If CV=View when the graph was constructed, then this method will be effective if the indices of the graph were already contiguous. In this case, the indices are left untouched and internal storage for the graph is minimized.
Graph().StorageOptimized()==true, if successful.
| int Epetra_CrsMatrix::PutScalar | ( | double | ScalarConstant | ) |
Initialize all values in the matrix with constant value.
| ScalarConstant | - (In) Value to use. |
| const Epetra_Map& Epetra_CrsMatrix::RangeMap | ( | ) | const [inline] |
Returns the Epetra_Map object associated with the range of this matrix operator.
| int Epetra_CrsMatrix::ReplaceColMap | ( | const Epetra_BlockMap & | newmap | ) |
Replaces the current ColMap with the user-specified map object.
Replaces the current ColMap with the user-specified map object, but only if currentmap->PointSameAs(newmap) is true. This is a collective function. Returns 0 if map is replaced, -1 if not.
| int Epetra_CrsMatrix::ReplaceDiagonalValues | ( | const Epetra_Vector & | 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_Vector will be ignored and the return code will be set to 1.
The Epetra_Map associated with the input Epetra_Vector must be compatible with the RowMap of the matrix.
| Diagonal | - (In) New values to be placed in the main diagonal. |
| virtual int Epetra_CrsMatrix::ReplaceGlobalValues | ( | int | GlobalRow, | |
| int | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) | [virtual] |
Replace specified existing values with this list of entries for a given global row of the matrix.
| GlobalRow | - (In) Row number (in global coordinates) to put elements. | |
| NumEntries | - (In) Number of entries. | |
| Values | - (In) Values to enter. | |
| Indices | - (In) Global column indices corresponding to values. |
Reimplemented in Epetra_FECrsMatrix.
| int Epetra_CrsMatrix::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.
| 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. |
Reimplemented in Epetra_OskiMatrix.
| int Epetra_CrsMatrix::ReplaceRowMap | ( | const Epetra_BlockMap & | newmap | ) |
Replaces the current RowMap with the user-specified map object.
Replaces the current RowMap with the user-specified map object, but only if currentmap->PointSameAs(newmap) is true. This is a collective function. Returns 0 if map is replaced, -1 if not.
| int Epetra_CrsMatrix::RightScale | ( | const Epetra_Vector & | x | ) | [virtual] |
Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x.
The this matrix will be scaled such that A(i,j) = x(j)*A(i,j) where i denotes the global row number of A and j denotes the global column number of A.
| x | - (In) The Epetra_Vector used for scaling this. |
Implements Epetra_RowMatrix.
| int Epetra_CrsMatrix::Scale | ( | double | ScalarConstant | ) |
Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A).
| ScalarConstant | - (In) Value to use. |
| int Epetra_CrsMatrix::SetUseTranspose | ( | bool | UseTranspose_in | ) | [inline, virtual] |
If set true, transpose of this operator will be applied.
This flag allows the transpose of the given operator to be used implicitly. Setting this flag affects only the Apply() and ApplyInverse() methods. If the implementation of this interface does not support transpose use, this method should return a value of -1.
| UseTranspose | - (In) If true, multiply by the transpose of operator, otherwise just use operator. |
Implements Epetra_Operator.
| int Epetra_CrsMatrix::Solve | ( | bool | Upper, | |
| bool | Trans, | |||
| bool | UnitDiagonal, | |||
| const Epetra_MultiVector & | X, | |||
| Epetra_MultiVector & | Y | |||
| ) | const [virtual] |
Returns the result of a local solve using the Epetra_CrsMatrix a Epetra_MultiVector X in Y.
This method solves a triangular system of equations asynchronously on each processor.
| Upper | - (In) If true, solve Uy = x, otherwise solve Ly = x. | |
| Trans | - (In) If true, solve transpose problem. | |
| UnitDiagonal | - (In) If true, assume diagonal is unit (whether it's stored or not). | |
| X | - (In) An Epetra_MultiVector of dimension NumVectors to solve for. | |
| Y | - (Out) An Epetra_MultiVector of dimension NumVectors containing result. |
Implements Epetra_RowMatrix.
Reimplemented in Epetra_OskiMatrix.
| int Epetra_CrsMatrix::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.
This method solves a triangular system of equations asynchronously on each processor.
| Upper | - (In) If true, solve Uy = x, otherwise solve Ly = x. | |
| Trans | - (In) If true, solve transpose problem. | |
| UnitDiagonal | - (In) If true, assume diagonal is unit (whether it's stored or not). | |
| x | - (In) An Epetra_Vector to solve for. | |
| y | - (Out) An Epetra_Vector containing result. |
Reimplemented in Epetra_OskiMatrix.
| int Epetra_CrsMatrix::SortGhostsAssociatedWithEachProcessor | ( | bool | Flag | ) | [inline] |
Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order.
To be compliant with AztecOO, FillComplete() already locally orders ghostnodes such that information received from processor k has a lower local numbering than information received from processor j if k is less than j. SortGhostsAssociatedWithEachProcessor(True) further forces FillComplete() to locally number all ghostnodes received from processor k in ascending order. That is, the local numbering of b is less than c if the global numbering of b is less than c and if both b and c are owned by the same processor. This is done to be compliant with some limited block features within ML. In particular, some ML features require that a block structure of the matrix be maintained even within the ghost variables. Always returns 0.
| virtual int Epetra_CrsMatrix::SumIntoGlobalValues | ( | int | GlobalRow, | |
| int | NumEntries, | |||
| double * | Values, | |||
| int * | Indices | |||
| ) | [virtual] |
Add this list of entries to existing values for a given global row of the matrix.
| GlobalRow | - (In) Row number (in global coordinates) to put elements. | |
| NumEntries | - (In) Number of entries. | |
| Values | - (In) Values to enter. | |
| Indices | - (In) Global column indices corresponding to values. |
Reimplemented in Epetra_FECrsMatrix.
| int Epetra_CrsMatrix::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.
| 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. |
Reimplemented in Epetra_OskiMatrix.
1.4.7