|
Epetra Package Browser (Single Doxygen Collection) Development
|
Epetra Finite-Element CrsMatrix. More...
#include <Epetra_FECrsMatrix.h>

Public Types | |
| enum | { ROW_MAJOR = 0 , COLUMN_MAJOR = 3 } |
Public Member Functions | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int *NumEntriesPerRow, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int NumEntriesPerRow, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, int *NumEntriesPerRow, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, int NumEntriesPerRow, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_CrsGraph &Graph, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (Epetra_DataAccess CV, const Epetra_FECrsGraph &Graph, bool ignoreNonLocalEntries=false) | |
| Constructor. | |
| Epetra_FECrsMatrix (const Epetra_FECrsMatrix &src) | |
| Copy Constructor. | |
| virtual | ~Epetra_FECrsMatrix () |
| Destructor. | |
| Epetra_FECrsMatrix & | operator= (const Epetra_FECrsMatrix &src) |
| Assignment operator. | |
| void | Print (std::ostream &os) const |
| Print method. | |
| int | SumIntoGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
| override base-class Epetra_CrsMatrix::SumIntoGlobalValues method | |
| int | SumIntoGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
| int | InsertGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
| override base-class Epetra_CrsMatrix::InsertGlobalValues method | |
| int | InsertGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
| int | InsertGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
| override base-class Epetra_CrsMatrix::InsertGlobalValues method | |
| int | InsertGlobalValues (long long GlobalRow, int NumEntries, double *Values, long long *Indices) |
| int | ReplaceGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
| override base-class Epetra_CrsMatrix::ReplaceGlobalValues method | |
| int | ReplaceGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
| int | SumIntoGlobalValues (int numIndices, const int *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Sum a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations. | |
| int | SumIntoGlobalValues (int numIndices, const long long *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | SumIntoGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Sum a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations. | |
| int | SumIntoGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | SumIntoGlobalValues (int numIndices, const int *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Sum C-style table (double-pointer, or list of lists) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations. | |
| int | SumIntoGlobalValues (int numIndices, const long long *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | SumIntoGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Sum C-style table (double-pointer, or list of lists) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations. | |
| int | SumIntoGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | InsertGlobalValues (int numIndices, const int *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Insert a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix. | |
| int | InsertGlobalValues (int numIndices, const long long *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | InsertGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Insert a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix. | |
| int | InsertGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | InsertGlobalValues (int numIndices, const int *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Insert a C-style table (double-pointer, or list of lists) of coefficients into the matrix. | |
| int | InsertGlobalValues (int numIndices, const long long *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | InsertGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Insert a C-style table (double-pointer, or list of lists) of coefficients into the matrix. | |
| int | InsertGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | ReplaceGlobalValues (int numIndices, const int *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Copy a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations. | |
| int | ReplaceGlobalValues (int numIndices, const long long *indices, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | ReplaceGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Copy Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations. | |
| int | ReplaceGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | ReplaceGlobalValues (int numIndices, const int *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Copy C-style table (double-pointer, or list of lists) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations. | |
| int | ReplaceGlobalValues (int numIndices, const long long *indices, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | ReplaceGlobalValues (int numRows, const int *rows, int numCols, const int *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| Copy C-style table (double-pointer, or list of lists) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations. | |
| int | ReplaceGlobalValues (int numRows, const long long *rows, int numCols, const long long *cols, const double *const *values, int format=Epetra_FECrsMatrix::ROW_MAJOR) |
| int | SumIntoGlobalValues (const Epetra_IntSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Sum a square structurally-symmetric sub-matrix into the global matrix. | |
| int | SumIntoGlobalValues (const Epetra_LongLongSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | SumIntoGlobalValues (const Epetra_IntSerialDenseVector &rows, const Epetra_IntSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Sum a general sub-matrix into the global matrix. | |
| int | SumIntoGlobalValues (const Epetra_LongLongSerialDenseVector &rows, const Epetra_LongLongSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | InsertGlobalValues (const Epetra_IntSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Insert a square structurally-symmetric sub-matrix into the global matrix. | |
| int | InsertGlobalValues (const Epetra_LongLongSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | InsertGlobalValues (const Epetra_IntSerialDenseVector &rows, const Epetra_IntSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Insert a general sub-matrix into the global matrix. | |
| int | InsertGlobalValues (const Epetra_LongLongSerialDenseVector &rows, const Epetra_LongLongSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | ReplaceGlobalValues (const Epetra_IntSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Use a square structurally-symmetric sub-matrix to replace existing values in the global matrix. | |
| int | ReplaceGlobalValues (const Epetra_LongLongSerialDenseVector &indices, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | ReplaceGlobalValues (const Epetra_IntSerialDenseVector &rows, const Epetra_IntSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| Use a general sub-matrix to replace existing values. | |
| int | ReplaceGlobalValues (const Epetra_LongLongSerialDenseVector &rows, const Epetra_LongLongSerialDenseVector &cols, const Epetra_SerialDenseMatrix &values, int format=Epetra_FECrsMatrix::COLUMN_MAJOR) |
| int | GlobalAssemble (bool callFillComplete=true, Epetra_CombineMode combineMode=Add, bool save_off_and_reuse_map_exporter=false) |
| Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time. | |
| int | GlobalAssemble (const Epetra_Map &domain_map, const Epetra_Map &range_map, bool callFillComplete=true, Epetra_CombineMode combineMode=Add, bool save_off_and_reuse_map_exporter=false) |
| Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time. | |
| void | setIgnoreNonLocalEntries (bool flag) |
| Set whether or not non-local data values should be ignored. | |
| 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. | |
| 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. | |
| 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 &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(). | |
| 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(). | |
| 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(). | |
| 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(). | |
| Epetra_CrsMatrix (const Epetra_CrsMatrix &Matrix) | |
| Copy constructor. | |
| virtual | ~Epetra_CrsMatrix () |
| Epetra_CrsMatrix Destructor. | |
| 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). | |
| 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. | |
| 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, const double *Values, const int *Indices) |
| Replace current values with this list of entries for a given local row of the matrix. | |
| 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. | |
| int | ReplaceDiagonalValues (const Epetra_Vector &Diagonal) |
| Replaces diagonal values of the matrix with those in the user-provided vector. | |
| 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. | |
| 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 | 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. | |
| 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 | 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. | |
| 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 | 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. | |
| int | ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values) const |
| Returns a view of the specified global row values via pointers to internal data. | |
| 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. | |
| 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. | |
| 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. | |
| int | NumMyDiagonals () const |
| Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. | |
| 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. | |
| 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. | |
| long long | IndexBase64 () const |
| 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. | |
| int | ReplaceDomainMapAndImporter (const Epetra_Map &NewDomainMap, const Epetra_Import *NewImporter) |
| Replaces the current DomainMap & Importer with the user-specified map object. | |
| int | RemoveEmptyProcessesInPlace (const Epetra_BlockMap *NewMap) |
| Remove processes owning zero rows from the Maps and their communicator. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| 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 |
| 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. | |
| int | NumMyRowEntries (int MyRow, int &NumEntries) const |
| Return the current number of values stored for the specified local row. | |
| const Epetra_BlockMap & | Map () const |
| Map() method inherited from Epetra_DistObject. | |
| 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. | |
| double * | operator[] (int Loc) |
| Inlined bracket operator for fast access to data. (Const and Non-const versions) | |
| double * | operator[] (int Loc) const |
| int | ExtractCrsDataPointers (int *&IndexOffset, int *&Indices, double *&Values_in) const |
| Returns internal data pointers associated with Crs matrix format. | |
| Epetra_IntSerialDenseVector & | ExpertExtractIndexOffset () |
| Returns a reference to the Epetra_IntSerialDenseVector used to hold the local IndexOffsets (CRS rowptr) | |
| Epetra_IntSerialDenseVector & | ExpertExtractIndices () |
| Returns a reference to the Epetra_IntSerialDenseVector used to hold the local All_Indices (CRS colind) | |
| double *& | ExpertExtractValues () |
| Returns a reference to the double* used to hold the values array. | |
| 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. | |
| int | ExpertMakeUniqueCrsGraphData () |
| Makes sure this matrix has a unique CrsGraphData object. | |
| int | SortGhostsAssociatedWithEachProcessor (bool Flag) |
| Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order. | |
| Epetra_DistObject (const Epetra_BlockMap &Map) | |
| Basic Epetra_DistObject constuctor. | |
| 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. | |
| 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. | |
| 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. | |
| 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. | |
| const Epetra_BlockMap & | Map () const |
| Returns the address of the Epetra_BlockMap for this multi-vector. | |
| const Epetra_Comm & | Comm () 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. | |
| virtual int | ReportError (const std::string Message, int ErrorCode) const |
| Error reporting method. | |
| Epetra_Object (int TracebackModeIn=-1, bool set_label=true) | |
| Epetra_Object Constructor. | |
| Epetra_Object (const char *const Label, int TracebackModeIn=-1) | |
| Epetra_Object Constructor. | |
| Epetra_Object (const Epetra_Object &Object) | |
| Epetra_Object Copy Constructor. | |
| virtual | ~Epetra_Object () |
| Epetra_Object Destructor. | |
| virtual void | SetLabel (const char *const Label) |
| Epetra_Object Label definition using char *. | |
| virtual | ~Epetra_SrcDistObject () |
| Epetra_SrcDistObject destructor. | |
| 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. | |
| Epetra_CompObject & | operator= (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_Flops * | GetFlopCounter () 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. | |
| Epetra_BLAS (void) | |
| Epetra_BLAS Constructor. | |
| Epetra_BLAS (const Epetra_BLAS &BLAS) | |
| Epetra_BLAS Copy Constructor. | |
| 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) | |
| virtual | ~Epetra_RowMatrix () |
| Destructor. | |
| virtual | ~Epetra_Operator () |
| Destructor. | |
Private Types | |
| enum | { SUMINTO = 0 , REPLACE = 1 , INSERT = 2 } |
Private Member Functions | |
| void | DeleteMemory () |
| template<typename int_type> | |
| int | InputGlobalValues (int numRows, const int_type *rows, int numCols, const int_type *cols, const double *const *values, int format, int mode) |
| template<typename int_type> | |
| int | InputGlobalValues (int numRows, const int_type *rows, int numCols, const int_type *cols, const double *values, int format, int mode) |
| template<typename int_type> | |
| int | InputNonlocalGlobalValues (int_type row, int numCols, const int_type *cols, const double *values, int mode) |
| template<typename int_type> | |
| int | InputGlobalValues_RowMajor (int numRows, const int_type *rows, int numCols, const int_type *cols, const double *values, int mode) |
| template<typename int_type> | |
| int | InsertNonlocalRow (int_type row, typename std::vector< int_type >::iterator offset) |
| template<typename int_type> | |
| int | InputNonlocalValue (int rowoffset, int_type col, double value, int mode) |
| template<typename int_type> | |
| std::vector< int_type > & | nonlocalRows () |
| template<typename int_type> | |
| std::vector< std::vector< int_type > > & | nonlocalCols () |
| template<typename int_type> | |
| int | SumIntoGlobalValues (int_type GlobalRow, int NumEntries, const double *values, const int_type *Indices) |
| template<typename int_type> | |
| int | GlobalAssemble (const Epetra_Map &domain_map, const Epetra_Map &range_map, bool callFillComplete=true, Epetra_CombineMode combineMode=Add, bool save_off_and_reuse_map_exporter=false) |
| template<> | |
| std::vector< int > & | nonlocalRows () |
| template<> | |
| std::vector< long long > & | nonlocalRows () |
| template<> | |
| std::vector< std::vector< int > > & | nonlocalCols () |
| template<> | |
| std::vector< std::vector< long long > > & | nonlocalCols () |
Private Attributes | |
| long long | myFirstRow_ |
| int | myNumRows_ |
| bool | ignoreNonLocalEntries_ |
| std::vector< int > | nonlocalRows_int_ |
| std::vector< std::vector< int > > | nonlocalCols_int_ |
| std::vector< long long > | nonlocalRows_LL_ |
| std::vector< std::vector< long long > > | nonlocalCols_LL_ |
| std::vector< std::vector< double > > | nonlocalCoefs_ |
| std::vector< double > | workData_ |
| std::vector< const double * > | workData2d_ |
| int | workDataLength_ |
| bool | useNonlocalMatrix_ |
| Epetra_CrsMatrix * | nonlocalMatrix_ |
| Epetra_Map * | sourceMap_ |
| Epetra_Map * | colMap_ |
| Epetra_Export * | exporter_ |
| Epetra_CrsMatrix * | tempMat_ |
Additional Inherited Members | |
| static void | SetTracebackMode (int TracebackModeValue) |
| Set the value of the Epetra_Object error traceback report mode. | |
| 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 int | TracebackMode |
| 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 () |
| 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. | |
| std::string | toString (const int &x) const |
| std::string | toString (const long long &x) const |
| std::string | toString (const double &x) const |
| 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_MultiVector * | ImportVector_ |
| Epetra_MultiVector * | ExportVector_ |
| Epetra_DataAccess | CV_ |
| bool | squareFillCompleteCalled_ |
| Epetra_BlockMap | Map_ |
| const Epetra_Comm * | Comm_ |
| char * | Exports_ |
| char * | Imports_ |
| int | LenExports_ |
| int | LenImports_ |
| int * | Sizes_ |
| Epetra_Flops * | FlopCounter_ |
Epetra Finite-Element CrsMatrix.
This class provides the ability to input finite-element style sub-matrix data, including sub-matrices with non-local rows (which could correspond to shared finite-element nodes for example). This class inherits Epetra_CrsMatrix, and so all Epetra_CrsMatrix functionality is also available.
It is intended that this class will be used as follows:
Sub-matrix data, which is assumed to be a rectangular 'table' of coefficients accompanied by 'scatter-indices', can be provided in three forms:
In all cases, a "format" parameter specifies whether the data is laid out in row-major or column-major order (i.e., whether coefficients for a row lie contiguously or whether coefficients for a column lie contiguously). See the documentation for the methods SumIntoGlobalValues() and ReplaceGlobalValues().
Important notes:
Definition at line 120 of file Epetra_FECrsMatrix.h.
| anonymous enum |
| Enumerator | |
|---|---|
| ROW_MAJOR | |
| COLUMN_MAJOR | |
Definition at line 167 of file Epetra_FECrsMatrix.h.
|
private |
| Enumerator | |
|---|---|
| SUMINTO | |
| REPLACE | |
| INSERT | |
Definition at line 748 of file Epetra_FECrsMatrix.h.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_Map & | RowMap, | ||
| int * | NumEntriesPerRow, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 54 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_Map & | RowMap, | ||
| int | NumEntriesPerRow, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 84 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_Map & | RowMap, | ||
| const Epetra_Map & | ColMap, | ||
| int * | NumEntriesPerRow, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 114 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_Map & | RowMap, | ||
| const Epetra_Map & | ColMap, | ||
| int | NumEntriesPerRow, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 145 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_CrsGraph & | Graph, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 176 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | Epetra_DataAccess | CV, |
| const Epetra_FECrsGraph & | Graph, | ||
| bool | ignoreNonLocalEntries = false ) |
Constructor.
Definition at line 205 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix::Epetra_FECrsMatrix | ( | const Epetra_FECrsMatrix & | src | ) |
Copy Constructor.
Definition at line 235 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Destructor.
Definition at line 313 of file Epetra_FECrsMatrix.cpp.
| Epetra_FECrsMatrix & Epetra_FECrsMatrix::operator= | ( | const Epetra_FECrsMatrix & | src | ) |
Assignment operator.
Definition at line 260 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Print method.
Reimplemented from Epetra_CrsMatrix.
Definition at line 318 of file Epetra_FECrsMatrix.cpp.
|
virtual |
override base-class Epetra_CrsMatrix::SumIntoGlobalValues method
Reimplemented from Epetra_CrsMatrix.
Definition at line 798 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Reimplemented from Epetra_CrsMatrix.
Definition at line 810 of file Epetra_FECrsMatrix.cpp.
|
virtual |
override base-class Epetra_CrsMatrix::InsertGlobalValues method
Reimplemented from Epetra_CrsMatrix.
Definition at line 821 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Reimplemented from Epetra_CrsMatrix.
Definition at line 831 of file Epetra_FECrsMatrix.cpp.
|
virtual |
override base-class Epetra_CrsMatrix::InsertGlobalValues method
Reimplemented from Epetra_CrsMatrix.
Definition at line 841 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Reimplemented from Epetra_CrsMatrix.
Definition at line 850 of file Epetra_FECrsMatrix.cpp.
|
virtual |
override base-class Epetra_CrsMatrix::ReplaceGlobalValues method
Reimplemented from Epetra_CrsMatrix.
Definition at line 860 of file Epetra_FECrsMatrix.cpp.
|
virtual |
Reimplemented from Epetra_CrsMatrix.
Definition at line 869 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Sum a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | List, length numIndices*numIndices. Square sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 471 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 481 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Sum a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List of column-numbers (scatter-indices) for the sub-matrix. |
| values | List, length numRows*numCols. Rectangular sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 492 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 504 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Sum C-style table (double-pointer, or list of lists) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | Square sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 427 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 437 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Sum C-style table (double-pointer, or list of lists) of coefficients into the matrix, adding them to any coefficients that may already exist at the specified row/column locations.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List of column-numbers (scatter-indices) for the sub-matrix. |
| values | Rectangular sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 448 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 459 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Insert a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | List, length numIndices*numIndices. Square sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 738 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 749 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Insert a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List of column-numbers (scatter-indices) for the sub-matrix. |
| values | List, length numRows*numCols. Rectangular sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 760 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 772 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Insert a C-style table (double-pointer, or list of lists) of coefficients into the matrix.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | Square sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 692 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 703 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Insert a C-style table (double-pointer, or list of lists) of coefficients into the matrix.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List of column-numbers (scatter-indices) for the sub-matrix. |
| values | Rectangular sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 714 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 726 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Copy a Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | List, length numIndices*numIndices. Square sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 923 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 933 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Copy Fortran-style table (single-dimensional packed-list) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List, of column-numbers (scatter-indices) for the sub-matrix. |
| values | List, length numRows*numCols. Rectangular sub-matrix of coefficients, packed in a 1-D array. Data is packed either contiguously by row or by column, specified by the final parameter 'format'. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is COLUMN_MAJOR. |
Definition at line 944 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double * | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 955 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numIndices, |
| const int * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Copy C-style table (double-pointer, or list of lists) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations.
| numIndices | Number of rows (and columns) in the sub-matrix. |
| indices | List of scatter-indices (rows and columns) for the sub-matrix. |
| values | Square sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 879 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numIndices, |
| const long long * | indices, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 889 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numRows, |
| const int * | rows, | ||
| int | numCols, | ||
| const int * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Copy C-style table (double-pointer, or list of lists) of coefficients into the matrix, replacing any coefficients that may already exist at the specified row/column locations.
| numRows | Number of rows in the sub-matrix. |
| rows | List of row-numbers (scatter-indices) for the sub-matrix. |
| numCols | Number of columns in the sub-matrix. |
| cols | List of column-numbers (scatter-indices) for the sub-matrix. |
| values | Rectangular sub-matrix of coefficients, provided in a 2-D array, or double-pointer. |
| format | Specifies whether the data in 'values' is packed in column-major or row-major order. Valid values are Epetra_FECrsMatrix::ROW_MAJOR or Epetra_FECrsMatrix::COLUMN_MAJOR. This is an optional parameter, default value is ROW_MAJOR. |
Definition at line 900 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | int | numRows, |
| const long long * | rows, | ||
| int | numCols, | ||
| const long long * | cols, | ||
| const double *const * | values, | ||
| int | format = Epetra_FECrsMatrix::ROW_MAJOR ) |
Definition at line 911 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | const Epetra_IntSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Sum a square structurally-symmetric sub-matrix into the global matrix.
For non-square sub-matrices, see the other overloading of this method.
| indices | List of scatter-indices. indices.Length() must be the same as values.M() and values.N(). |
| values | Sub-matrix of coefficients. Must be square. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 516 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 530 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | const Epetra_IntSerialDenseVector & | rows, |
| const Epetra_IntSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Sum a general sub-matrix into the global matrix.
For square structurally-symmetric sub-matrices, see the other overloading of this method.
| rows | List of row-indices. rows.Length() must be the same as values.M(). |
| cols | List of column-indices. cols.Length() must be the same as values.N(). |
| values | Sub-matrix of coefficients. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 598 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::SumIntoGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | rows, |
| const Epetra_LongLongSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 613 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | const Epetra_IntSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Insert a square structurally-symmetric sub-matrix into the global matrix.
For non-square sub-matrices, see the other overloading of this method.
| indices | List of scatter-indices. indices.Length() must be the same as values.M() and values.N(). |
| values | Sub-matrix of coefficients. Must be square. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 544 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 557 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | const Epetra_IntSerialDenseVector & | rows, |
| const Epetra_IntSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Insert a general sub-matrix into the global matrix.
For square structurally-symmetric sub-matrices, see the other overloading of this method.
| rows | List of row-indices. rows.Length() must be the same as values.M(). |
| cols | List of column-indices. cols.Length() must be the same as values.N(). |
| values | Sub-matrix of coefficients. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 629 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::InsertGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | rows, |
| const Epetra_LongLongSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 644 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | const Epetra_IntSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Use a square structurally-symmetric sub-matrix to replace existing values in the global matrix.
For non-square sub-matrices, see the other overloading of this method.
| indices | List of scatter-indices. indices.Length() must be the same as values.M() and values.N(). |
| values | Sub-matrix of coefficients. Must be square. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 571 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | indices, |
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 584 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | const Epetra_IntSerialDenseVector & | rows, |
| const Epetra_IntSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Use a general sub-matrix to replace existing values.
For square structurally-symmetric sub-matrices, see the other overloading of this method.
| rows | List of row-indices. rows.Length() must be the same as values.M(). |
| cols | List of column-indices. cols.Length() must be the same as values.N(). |
| values | Sub-matrix of coefficients. |
| format | Optional format specifier, defaults to COLUMN_MAJOR. |
Definition at line 660 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::ReplaceGlobalValues | ( | const Epetra_LongLongSerialDenseVector & | rows, |
| const Epetra_LongLongSerialDenseVector & | cols, | ||
| const Epetra_SerialDenseMatrix & | values, | ||
| int | format = Epetra_FECrsMatrix::COLUMN_MAJOR ) |
Definition at line 676 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::GlobalAssemble | ( | bool | callFillComplete = true, |
| Epetra_CombineMode | combineMode = Add, | ||
| bool | save_off_and_reuse_map_exporter = false ) |
Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time.
Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method – every processor must enter it before any will complete it.
NOTE***: When GlobalAssemble() calls FillComplete(), it passes the arguments 'DomainMap()' and 'RangeMap()', which are the map attributes held by the base-class CrsMatrix and its graph. If a rectangular matrix is being assembled, the domain-map and range-map must be specified by calling the other overloading of this method. Otherwise, GlobalAssemble() has no way of knowing what these maps should really be.
| callFillComplete | option argument, defaults to true. Determines whether GlobalAssemble() internally calls the FillComplete() method on this matrix. |
Definition at line 966 of file Epetra_FECrsMatrix.cpp.
| int Epetra_FECrsMatrix::GlobalAssemble | ( | const Epetra_Map & | domain_map, |
| const Epetra_Map & | range_map, | ||
| bool | callFillComplete = true, | ||
| Epetra_CombineMode | combineMode = Add, | ||
| bool | save_off_and_reuse_map_exporter = false ) |
Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time.
Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method – every processor must enter it before any will complete it.
NOTE***: When GlobalAssemble() (the other overloading of this method) calls FillComplete(), it passes the arguments 'DomainMap()' and 'RangeMap()', which are the map attributes already held by the base-class CrsMatrix and its graph. If a rectangular matrix is being assembled, the domain-map and range-map must be specified. Otherwise, GlobalAssemble() has no way of knowing what these maps should really be.
| domain_map | user-supplied domain map for this matrix |
| range_map | user-supplied range map for this matrix |
| callFillComplete | option argument, defaults to true. Determines whether GlobalAssemble() internally calls the FillComplete() method on this matrix. |
Definition at line 1123 of file Epetra_FECrsMatrix.cpp.
|
inline |
Set whether or not non-local data values should be ignored.
By default, non-local data values are NOT ignored.
Definition at line 741 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 393 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1218 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1268 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1299 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1154 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1366 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 1386 of file Epetra_FECrsMatrix.cpp.
|
private |
|
private |
|
private |
Definition at line 784 of file Epetra_FECrsMatrix.cpp.
|
private |
Definition at line 974 of file Epetra_FECrsMatrix.cpp.
|
inlineprivate |
Definition at line 120 of file Epetra_FECrsMatrix.h.
|
inlineprivate |
Definition at line 120 of file Epetra_FECrsMatrix.h.
|
inlineprivate |
Definition at line 120 of file Epetra_FECrsMatrix.h.
|
inlineprivate |
Definition at line 120 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 785 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 786 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 788 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 791 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 792 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 795 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 796 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 802 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 806 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 807 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 808 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 810 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 811 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 813 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 814 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 815 of file Epetra_FECrsMatrix.h.
|
private |
Definition at line 816 of file Epetra_FECrsMatrix.h.