Epetra Package Browser (Single Doxygen Collection)
Development
|
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices. More...
#include <Epetra_CrsMatrix.h>
Public Member Functions | |
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) |
Public Member Functions inherited from Epetra_DistObject | |
Epetra_DistObject (const Epetra_BlockMap &Map) | |
Basic Epetra_DistObject constuctor. More... | |
Epetra_DistObject (const Epetra_BlockMap &Map, const char *const Label) | |
Epetra_DistObject (const Epetra_DistObject &Source) | |
Epetra_DistObject copy constructor. More... | |
virtual | ~Epetra_DistObject () |
Epetra_DistObject destructor. More... | |
int | Import (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Imports an Epetra_DistObject using the Epetra_Import object. More... | |
int | Import (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Imports an Epetra_DistObject using the Epetra_Export object. More... | |
int | Export (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Exports an Epetra_DistObject using the Epetra_Import object. More... | |
int | Export (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Exports an Epetra_DistObject using the Epetra_Export object. More... | |
const Epetra_BlockMap & | Map () const |
Returns the address of the Epetra_BlockMap for this multi-vector. More... | |
const Epetra_Comm & | Comm () const |
Returns the address of the Epetra_Comm for this multi-vector. More... | |
bool | DistributedGlobal () const |
Returns true if this multi-vector is distributed global, i.e., not local replicated. More... | |
Public Member Functions inherited from Epetra_Object | |
Epetra_Object (int TracebackModeIn=-1, bool set_label=true) | |
Epetra_Object Constructor. More... | |
Epetra_Object (const char *const Label, int TracebackModeIn=-1) | |
Epetra_Object Constructor. More... | |
Epetra_Object (const Epetra_Object &Object) | |
Epetra_Object Copy Constructor. More... | |
virtual | ~Epetra_Object () |
Epetra_Object Destructor. More... | |
virtual void | SetLabel (const char *const Label) |
Epetra_Object Label definition using char *. More... | |
virtual int | ReportError (const std::string Message, int ErrorCode) const |
Error reporting method. More... | |
Public Member Functions inherited from Epetra_SrcDistObject | |
virtual | ~Epetra_SrcDistObject () |
Epetra_SrcDistObject destructor. More... | |
Public Member Functions inherited from Epetra_CompObject | |
Epetra_CompObject & | operator= (const Epetra_CompObject &src) |
Epetra_CompObject () | |
Basic Epetra_CompObject constuctor. More... | |
Epetra_CompObject (const Epetra_CompObject &Source) | |
Epetra_CompObject copy constructor. More... | |
virtual | ~Epetra_CompObject () |
Epetra_CompObject destructor. More... | |
void | SetFlopCounter (const Epetra_Flops &FlopCounter_in) |
Set the internal Epetra_Flops() pointer. More... | |
void | SetFlopCounter (const Epetra_CompObject &CompObject) |
Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. More... | |
void | UnsetFlopCounter () |
Set the internal Epetra_Flops() pointer to 0 (no flops counted). More... | |
Epetra_Flops * | GetFlopCounter () const |
Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none. More... | |
void | ResetFlops () const |
Resets the number of floating point operations to zero for this multi-vector. More... | |
double | Flops () const |
Returns the number of floating point operations with this multi-vector. More... | |
void | UpdateFlops (int Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (long int Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (long long Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (double Flops_in) const |
Increment Flop count for this object. More... | |
void | UpdateFlops (float Flops_in) const |
Increment Flop count for this object. More... | |
Public Member Functions inherited from Epetra_BLAS | |
Epetra_BLAS (void) | |
Epetra_BLAS Constructor. More... | |
Epetra_BLAS (const Epetra_BLAS &BLAS) | |
Epetra_BLAS Copy Constructor. More... | |
virtual | ~Epetra_BLAS (void) |
Epetra_BLAS Destructor. More... | |
float | ASUM (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS one norm function (SASUM). More... | |
double | ASUM (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS one norm function (DASUM). More... | |
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). More... | |
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). More... | |
float | NRM2 (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS norm function (SNRM2). More... | |
double | NRM2 (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS norm function (DNRM2). More... | |
void | SCAL (const int N, const float ALPHA, float *X, const int INCX=1) const |
Epetra_BLAS vector scale function (SSCAL) More... | |
void | SCAL (const int N, const double ALPHA, double *X, const int INCX=1) const |
Epetra_BLAS vector scale function (DSCAL) More... | |
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) More... | |
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) More... | |
int | IAMAX (const int N, const float *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (ISAMAX) More... | |
int | IAMAX (const int N, const double *X, const int INCX=1) const |
Epetra_BLAS arg maximum of absolute value function (IDAMAX) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
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) More... | |
Public Member Functions inherited from Epetra_RowMatrix | |
virtual | ~Epetra_RowMatrix () |
Destructor. More... | |
Public Member Functions inherited from Epetra_Operator | |
virtual | ~Epetra_Operator () |
Destructor. More... | |
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 | 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. More... | |
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. More... | |
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(). More... | |
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(). More... | |
int | SortEntries () |
Sort column entries, row-by-row, in ascending order. More... | |
bool | Sorted () const |
If SortEntries() has been called, this query returns true, otherwise it returns false. More... | |
int | MergeRedundantEntries () |
Add entries that have the same column index. Remove redundant entries from list. More... | |
bool | NoRedundancies () const |
If MergeRedundantEntries() has been called, this query returns true, otherwise it returns false. More... | |
void | DeleteMemory () |
Protected Member Functions inherited from Epetra_DistObject | |
virtual int | DoTransfer (const Epetra_SrcDistObject &A, Epetra_CombineMode CombineMode, int NumSameIDs, int NumPermuteIDs, int NumRemoteIDs, int NumExportIDs, int *PermuteToLIDs, int *PermuteFromLIDs, int *RemoteLIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &LenImports, char *&Imports, Epetra_Distributor &Distor, bool DoReverse, const Epetra_OffsetIndex *Indexor) |
Perform actual transfer (redistribution) of data across memory images, using Epetra_Distributor object. More... | |
Protected Member Functions inherited from Epetra_Object | |
std::string | toString (const int &x) const |
std::string | toString (const long long &x) const |
std::string | toString (const double &x) const |
Protected Attributes | |
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_ |
Protected Attributes inherited from Epetra_DistObject | |
Epetra_BlockMap | Map_ |
const Epetra_Comm * | Comm_ |
char * | Exports_ |
char * | Imports_ |
int | LenExports_ |
int | LenImports_ |
int * | Sizes_ |
Protected Attributes inherited from Epetra_CompObject | |
Epetra_Flops * | FlopCounter_ |
Private Member Functions | |
int | Solve1 (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const |
int | Solve1 (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
template<typename int_type > | |
int | TInsertGlobalValues (int_type Row, int NumEntries, const double *values, const int_type *Indices) |
template<typename int_type > | |
int | TInsertGlobalValues (int_type Row, int NumEntries, double *values, int_type *Indices) |
template<typename int_type > | |
int | InsertValues (int Row, int NumEntries, const double *values, const int_type *Indices) |
template<typename int_type > | |
int | InsertValues (int Row, int NumEntries, double *values, int_type *Indices) |
template<typename int_type > | |
int | TReplaceGlobalValues (int_type Row, int NumEntries, const double *srcValues, const int_type *Indices) |
template<typename int_type > | |
int | TSumIntoGlobalValues (int_type Row, int NumEntries, const double *srcValues, const int_type *Indices) |
template<typename int_type > | |
int | ExtractGlobalRowCopy (int_type Row, int Length, int &NumEntries, double *values, int_type *Indices) const |
template<typename int_type > | |
int | ExtractGlobalRowCopy (int_type Row, int Length, int &NumEntries, double *values) const |
template<typename int_type > | |
int | ExtractGlobalRowView (int_type Row, int &NumEntries, double *&values, int_type *&Indices) const |
template<typename int_type > | |
int | ExtractGlobalRowView (int_type Row, int &NumEntries, double *&values) const |
template<typename int_type > | |
int | TCopyAndPermuteCrsMatrix (const Epetra_CrsMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
template<typename int_type > | |
int | TCopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
template<typename int_type > | |
int | TUnpackAndCombine (const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor) |
template<class TransferType > | |
void | FusedTransfer (const Epetra_CrsMatrix &SourceMatrix, const TransferType &RowTransfer, const TransferType *DomainTransfer, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator) |
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. More... | |
Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int NumEntriesPerRow, bool StaticProfile=false) | |
Epetra_CrsMatrix constructor with fixed number of indices per row. More... | |
Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, const int *NumEntriesPerRow, bool StaticProfile=false) | |
Epetra_CrsMatrix constructor with variable number of indices per row. More... | |
Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, int NumEntriesPerRow, bool StaticProfile=false) | |
Epetra_CrsMatrix constuctor with fixed number of indices per row. More... | |
Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_CrsGraph &Graph) | |
Construct a matrix using an existing Epetra_CrsGraph object. More... | |
Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false) | |
Epetra CrsMatrix constructor that also fuses Import and FillComplete(). More... | |
Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Import *DomainImporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator) | |
Epetra CrsMatrix constructor that also fuses Import and FillComplete(). More... | |
Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false) | |
Epetra CrsMatrix constructor that also fuses Ex[prt and FillComplete(). More... | |
Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Export *DomainExporter, const Epetra_Map *DomainMap, const Epetra_Map *RangeMap, bool RestrictCommunicator) | |
Epetra CrsMatrix constructor that also fuses Ex[prt and FillComplete(). More... | |
Epetra_CrsMatrix (const Epetra_CrsMatrix &Matrix) | |
Copy constructor. More... | |
virtual | ~Epetra_CrsMatrix () |
Epetra_CrsMatrix Destructor. More... | |
Insertion/Replace/SumInto methods | |
Epetra_CrsMatrix & | operator= (const Epetra_CrsMatrix &src) |
Assignment operator. More... | |
int | PutScalar (double ScalarConstant) |
Initialize all values in the matrix with constant value. More... | |
int | Scale (double ScalarConstant) |
Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A). More... | |
virtual int | InsertGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
Insert a list of elements in a given global row of the matrix. More... | |
virtual int | InsertGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
virtual int | InsertGlobalValues (int GlobalRow, int NumEntries, double *Values, int *Indices) |
Insert a list of elements in a given global row of the matrix. More... | |
virtual int | InsertGlobalValues (long long GlobalRow, int NumEntries, double *Values, long long *Indices) |
virtual int | ReplaceGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
Replace specified existing values with this list of entries for a given global row of the matrix. More... | |
virtual int | ReplaceGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
virtual int | SumIntoGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices) |
Add this list of entries to existing values for a given global row of the matrix. More... | |
virtual int | SumIntoGlobalValues (long long GlobalRow, int NumEntries, const double *Values, const long long *Indices) |
int | InsertMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices) |
Insert a list of elements in a given local row of the matrix. More... | |
int | InsertMyValues (int MyRow, int NumEntries, double *Values, int *Indices) |
Insert a list of elements in a given local row of the matrix. More... | |
int | ReplaceMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices) |
Replace current values with this list of entries for a given local row of the matrix. More... | |
int | SumIntoMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices) |
Add this list of entries to existing values for a given local row of the matrix. More... | |
int | ReplaceDiagonalValues (const Epetra_Vector &Diagonal) |
Replaces diagonal values of the matrix with those in the user-provided vector. More... | |
Transformation methods | |
int | FillComplete (bool OptimizeDataStorage=true) |
Signal that data entry is complete. Perform transformations to local index space. More... | |
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. More... | |
int | OptimizeStorage () |
Make consecutive row index sections contiguous, minimize internal storage used for constructing graph. More... | |
int | MakeDataContiguous () |
Eliminates memory that is used for construction. Make consecutive row index sections contiguous. More... | |
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. More... | |
int | ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values, long long *Indices) const |
int | ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const |
Returns a copy of the specified local row in user-provided arrays. More... | |
int | ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values) const |
Returns a copy of the specified global row values in user-provided array. More... | |
int | ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values) const |
int | ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values) const |
Returns a copy of the specified local row values in user-provided array. More... | |
int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const |
Returns a copy of the main diagonal in a user-provided vector. More... | |
int | ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values, int *&Indices) const |
Returns a view of the specified global row values via pointers to internal data. More... | |
int | ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values, long long *&Indices) const |
int | ExtractMyRowView (int MyRow, int &NumEntries, double *&Values, int *&Indices) const |
Returns a view of the specified local row values via pointers to internal data. More... | |
int | ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values) const |
Returns a view of the specified global row values via pointers to internal data. More... | |
int | ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values) const |
int | ExtractMyRowView (int MyRow, int &NumEntries, double *&Values) const |
Returns a view of the specified local row values via pointers to internal data. More... | |
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. More... | |
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. More... | |
int | Multiply1 (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
int | Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const |
Returns the result of a local solve using the Epetra_CrsMatrix on a Epetra_Vector x in y. More... | |
int | 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. More... | |
int | InvRowSums (Epetra_Vector &x) const |
Computes the inverse of the sum of absolute values of the rows of the Epetra_CrsMatrix, results returned in x. More... | |
int | InvRowMaxs (Epetra_Vector &x) const |
Computes the inverse of the max of absolute values of the rows of the Epetra_CrsMatrix, results returned in x. More... | |
int | LeftScale (const Epetra_Vector &x) |
Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x. More... | |
int | InvColSums (Epetra_Vector &x) const |
Computes the inverse of the sum of absolute values of the columns of the Epetra_CrsMatrix, results returned in x. More... | |
int | InvColMaxs (Epetra_Vector &x) const |
Computes the max of absolute values of the columns of the Epetra_CrsMatrix, results returned in x. More... | |
int | RightScale (const Epetra_Vector &x) |
Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x. More... | |
Matrix Properties Query Methods | |
bool | Filled () const |
If FillComplete() has been called, this query returns true, otherwise it returns false. More... | |
bool | StorageOptimized () const |
If OptimizeStorage() has been called, this query returns true, otherwise it returns false. More... | |
bool | IndicesAreGlobal () const |
If matrix indices has not been transformed to local, this query returns true, otherwise it returns false. More... | |
bool | IndicesAreLocal () const |
If matrix indices has been transformed to local, this query returns true, otherwise it returns false. More... | |
bool | IndicesAreContiguous () const |
If matrix indices are packed into single array (done in OptimizeStorage()) return true, otherwise false. More... | |
bool | LowerTriangular () const |
If matrix is lower triangular in local index space, this query returns true, otherwise it returns false. More... | |
bool | UpperTriangular () const |
If matrix is upper triangular in local index space, this query returns true, otherwise it returns false. More... | |
bool | NoDiagonal () const |
If matrix has no diagonal entries in global index space, this query returns true, otherwise it returns false. More... | |
Attribute access functions | |
double | NormInf () const |
Returns the infinity norm of the global matrix. More... | |
double | NormOne () const |
Returns the one norm of the global matrix. More... | |
double | NormFrobenius () const |
Returns the frobenius norm of the global matrix. More... | |
int | NumGlobalNonzeros () const |
Returns the number of nonzero entries in the global matrix. More... | |
long long | NumGlobalNonzeros64 () const |
int | NumGlobalRows () const |
Returns the number of global matrix rows. More... | |
long long | NumGlobalRows64 () const |
int | NumGlobalCols () const |
Returns the number of global matrix columns. More... | |
long long | NumGlobalCols64 () const |
int | NumGlobalDiagonals () const |
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. More... | |
long long | NumGlobalDiagonals64 () const |
int | NumMyNonzeros () const |
Returns the number of nonzero entries in the calling processor's portion of the matrix. More... | |
int | NumMyRows () const |
Returns the number of matrix rows owned by the calling processor. More... | |
int | NumMyCols () const |
Returns the number of entries in the set of column-indices that appear on this processor. More... | |
int | NumMyDiagonals () const |
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. More... | |
int | NumGlobalEntries (long long Row) const |
Returns the current number of nonzero entries in specified global row on this processor. More... | |
int | NumAllocatedGlobalEntries (int Row) const |
Returns the allocated number of nonzero entries in specified global row on this processor. More... | |
int | MaxNumEntries () const |
Returns the maximum number of nonzero entries across all rows on this processor. More... | |
int | GlobalMaxNumEntries () const |
Returns the maximum number of nonzero entries across all rows on all processors. More... | |
int | NumMyEntries (int Row) const |
Returns the current number of nonzero entries in specified local row on this processor. More... | |
int | NumAllocatedMyEntries (int Row) const |
Returns the allocated number of nonzero entries in specified local row on this processor. More... | |
int | IndexBase () const |
Returns the index base for row and column indices for this graph. More... | |
long long | IndexBase64 () const |
bool | StaticGraph () |
Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable. More... | |
const Epetra_CrsGraph & | Graph () const |
Returns a reference to the Epetra_CrsGraph object associated with this matrix. More... | |
const Epetra_Map & | RowMap () const |
Returns the Epetra_Map object associated with the rows of this matrix. More... | |
int | ReplaceRowMap (const Epetra_BlockMap &newmap) |
Replaces the current RowMap with the user-specified map object. More... | |
bool | HaveColMap () const |
Returns true if we have a well-defined ColMap, and returns false otherwise. More... | |
int | ReplaceColMap (const Epetra_BlockMap &newmap) |
Replaces the current ColMap with the user-specified map object. More... | |
int | ReplaceDomainMapAndImporter (const Epetra_Map &NewDomainMap, const Epetra_Import *NewImporter) |
Replaces the current DomainMap & Importer with the user-specified map object. More... | |
int | RemoveEmptyProcessesInPlace (const Epetra_BlockMap *NewMap) |
Remove processes owning zero rows from the Maps and their communicator. More... | |
const Epetra_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. More... | |
const Epetra_Map & | DomainMap () const |
Returns the Epetra_Map object associated with the domain of this matrix operator. More... | |
const Epetra_Map & | RangeMap () const |
Returns the Epetra_Map object associated with the range of this matrix operator. More... | |
const Epetra_Import * | Importer () const |
Returns the Epetra_Import object that contains the import operations for distributed operations. More... | |
const Epetra_Export * | Exporter () const |
Returns the Epetra_Export object that contains the export operations for distributed operations. More... | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this matrix. More... | |
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. More... | |
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. More... | |
long long | GRID64 (int LRID_in) const |
int | LCID (int GCID_in) const |
Returns the local column index for given global column index, returns -1 if no local column for this global column. More... | |
int | LCID (long long GCID_in) const |
int | GCID (int LCID_in) const |
Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column. More... | |
long long | GCID64 (int LCID_in) const |
bool | MyGRID (int GRID_in) const |
Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
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. More... | |
bool | MyGCID (int GCID_in) const |
Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
bool | MyGCID (long long GCID_in) const |
bool | MyLCID (int LCID_in) const |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
bool | MyGlobalRow (int GID) const |
Returns true of GID is owned by the calling processor, otherwise it returns false. More... | |
bool | MyGlobalRow (long long GID) const |
I/O Methods | |
virtual void | Print (std::ostream &os) const |
Print method. More... | |
Additional methods required to support the Epetra_Operator interface | |
const char * | Label () const |
Returns a character string describing the operator. More... | |
int | SetUseTranspose (bool UseTranspose_in) |
If set true, transpose of this operator will be applied. More... | |
int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More... | |
int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More... | |
bool | HasNormInf () const |
Returns true because this class can compute an Inf-norm. More... | |
bool | UseTranspose () const |
Returns the current UseTranspose setting. More... | |
const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this matrix operator. More... | |
const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this matrix operator. More... | |
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. More... | |
const Epetra_BlockMap & | Map () const |
Map() method inherited from Epetra_DistObject. More... | |
const Epetra_Map & | RowMatrixRowMap () const |
Returns the Epetra_Map object associated with the rows of this matrix. More... | |
const Epetra_Map & | RowMatrixColMap () const |
Returns the Epetra_Map object associated with columns of this matrix. More... | |
const Epetra_Import * | RowMatrixImporter () const |
Returns the Epetra_Import object that contains the import operations for distributed operations. More... | |
Inlined Operator Methods | |
double * | operator[] (int Loc) |
Inlined bracket operator for fast access to data. (Const and Non-const versions) More... | |
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. More... | |
Epetra_IntSerialDenseVector & | ExpertExtractIndexOffset () |
Returns a reference to the Epetra_IntSerialDenseVector used to hold the local IndexOffsets (CRS rowptr) More... | |
Epetra_IntSerialDenseVector & | ExpertExtractIndices () |
Returns a reference to the Epetra_IntSerialDenseVector used to hold the local All_Indices (CRS colind) More... | |
double *& | ExpertExtractValues () |
Returns a reference to the double* used to hold the values array. More... | |
int | ExpertStaticFillComplete (const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, const Epetra_Import *Importer=0, const Epetra_Export *Exporter=0, int NumMyDiagonals=-1) |
Performs a FillComplete on an object that aready has filled CRS data. More... | |
int | ExpertMakeUniqueCrsGraphData () |
Makes sure this matrix has a unique CrsGraphData object. More... | |
int | SortGhostsAssociatedWithEachProcessor (bool Flag) |
Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order. More... | |
Deprecated methods: These methods still work, but will be removed in a future version | |
const Epetra_Map & | ImportMap () const |
Use ColMap() instead. More... | |
int | TransformToLocal () |
Use FillComplete() instead. More... | |
int | TransformToLocal (const Epetra_Map *DomainMap, const Epetra_Map *RangeMap) |
Use FillComplete(const Epetra_Map& DomainMap, const Epetra_Map& RangeMap) instead. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from Epetra_Object | |
static void | SetTracebackMode (int TracebackModeValue) |
Set the value of the Epetra_Object error traceback report mode. More... | |
static int | GetTracebackMode () |
Get the value of the Epetra_Object error report mode. More... | |
static std::ostream & | GetTracebackStream () |
Get the output stream for error reporting. More... | |
Static Public Attributes inherited from Epetra_Object | |
static int | TracebackMode |
Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices.
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.
Definition at line 173 of file Epetra_CrsMatrix.h.
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. |
Definition at line 90 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 118 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 145 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 174 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 205 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix::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().
A common use case is to create an empty destination Epetra_CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call FillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.
Fusing redistribution and FillComplete() exposes potential optimizations. For example, it may make constructing the column map faster, and it may avoid intermediate unoptimized storage in the destination Epetra_CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.
The resulting matrix is fill complete (in the sense of Filled()) and has optimized storage (in the sense of StorageOptimized()). It the DomainMap is taken from the SourceMatrix, the RangeMap is presumed to be RowImporter.TargetMap() if not specified
SourceMatrix | [in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution. |
RowImporter | [in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the row Map of sourceMatrix. |
DomainMap | [in] The new domainMap for the new matrix. If not specified, then the DomainMap of the SourceMatrix is used. |
RangeMap | [in] The new rangeMap for the new matrix. If not specified, then RowImporter.TargetMap() is used. |
RestrictCommunicator | [in] Restricts the resulting communicator to active processes only. |
Definition at line 5127 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix::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().
A common use case is to create an empty destination Epetra_CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call FillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.
Fusing redistribution and FillComplete() exposes potential optimizations. For example, it may make constructing the column map faster, and it may avoid intermediate unoptimized storage in the destination Epetra_CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.
The resulting matrix is fill complete (in the sense of Filled()) and has optimized storage (in the sense of StorageOptimized()). It the DomainMap is taken from the SourceMatrix, the RangeMap is presumed to be RowImporter.TargetMap() if not specified
SourceMatrix | [in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution. |
RowImporter | [in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the row Map of sourceMatrix. |
DomainImporter | [in] The Import instance containing a precomputed redistribution plan (for the domain maps). The source Map of the Import must be the same as the domain Map of sourceMatrix. |
DomainMap | [in] The new domainMap for the new matrix. |
RangeMap | [in] The new rangeMap for the new matrix. |
RestrictCommunicator | [in] Restricts the resulting communicator to active processes only. |
Definition at line 5162 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix::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().
A common use case is to create an empty destination Epetra_CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call FillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.
Fusing redistribution and FillComplete() exposes potential optimizations. For example, it may make constructing the column map faster, and it may avoid intermediate unoptimized storage in the destination Epetra_CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.
The resulting matrix is fill complete (in the sense of Filled()) and has optimized storage (in the sense of StorageOptimized()). It the DomainMap is taken from the SourceMatrix, the RangeMap is presumed to be RowImporter.TargetMap() if not specified
SourceMatrix | [in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution. |
RowExporter | [in] The Export instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the row Map of sourceMatrix. |
DomainMap | [in] The new domainMap for the new matrix. If not specified, then the DomainMap of the SourceMatrix is used. |
RangeMap | [in] The new rangeMap for the new matrix. If not specified, then RowExporter.TargetMap() is used. |
RestrictCommunicator | [in] Restricts the resulting communicator to active processes only. |
Definition at line 5145 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix::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().
A common use case is to create an empty destination Epetra_CrsMatrix, redistribute from a source CrsMatrix (by an Import or Export operation), then call FillComplete() on the destination CrsMatrix. This constructor fuses these three cases, for an Import redistribution.
Fusing redistribution and FillComplete() exposes potential optimizations. For example, it may make constructing the column map faster, and it may avoid intermediate unoptimized storage in the destination Epetra_CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.
The resulting matrix is fill complete (in the sense of Filled()) and has optimized storage (in the sense of StorageOptimized()). It the DomainMap is taken from the SourceMatrix, the RangeMap is presumed to be RowImporter.TargetMap() if not specified
SourceMatrix | [in] The source matrix from which to import. The source of an Import must have a nonoverlapping distribution. |
RowExporter | [in] The Export instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the row Map of sourceMatrix. |
DomainExporter | [in] The Export instance containing a precomputed redistribution plan (for the domain map. The source Map of the Import must be the same as the domain Map of sourceMatrix. |
DomainMap | [in] The new domainMap for the new matrix. |
RangeMap | [in] The new rangeMap for the new matrix. |
RestrictCommunicator | [in] Restricts the resulting communicator to active processes only. |
Definition at line 5180 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix::Epetra_CrsMatrix | ( | const Epetra_CrsMatrix & | Matrix | ) |
Copy constructor.
Definition at line 234 of file Epetra_CrsMatrix.cpp.
|
virtual |
Epetra_CrsMatrix Destructor.
Definition at line 388 of file Epetra_CrsMatrix.cpp.
Epetra_CrsMatrix & Epetra_CrsMatrix::operator= | ( | const Epetra_CrsMatrix & | src | ) |
Assignment operator.
Definition at line 262 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::PutScalar | ( | double | ScalarConstant | ) |
Initialize all values in the matrix with constant value.
ScalarConstant | - (In) Value to use. |
Definition at line 487 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 504 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 540 of file Epetra_CrsMatrix.cpp.
|
virtual |
Reimplemented in Epetra_FECrsMatrix.
Definition at line 551 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 581 of file Epetra_CrsMatrix.cpp.
|
virtual |
Reimplemented in Epetra_FECrsMatrix.
Definition at line 592 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 852 of file Epetra_CrsMatrix.cpp.
|
virtual |
Reimplemented in Epetra_FECrsMatrix.
Definition at line 860 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1028 of file Epetra_CrsMatrix.cpp.
|
virtual |
Reimplemented in Epetra_FECrsMatrix.
Definition at line 1040 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::InsertMyValues | ( | int | MyRow, |
int | NumEntries, | ||
const double * | Values, | ||
const 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. |
Definition at line 604 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 625 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::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.
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. |
Definition at line 869 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::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.
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. |
Definition at line 1053 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1511 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::FillComplete | ( | bool | OptimizeDataStorage = true | ) |
Signal that data entry is complete. Perform transformations to local index space.
Definition at line 1140 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::FillComplete | ( | const Epetra_Map & | DomainMap, |
const Epetra_Map & | RangeMap, | ||
bool | OptimizeDataStorage = true |
||
) |
Signal that data entry is complete. Perform transformations to local index space.
Definition at line 1147 of file Epetra_CrsMatrix.cpp.
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.
Definition at line 1279 of file Epetra_CrsMatrix.cpp.
|
inline |
Eliminates memory that is used for construction. Make consecutive row index sections contiguous.
Definition at line 718 of file Epetra_CrsMatrix.h.
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.
\param GlobalRow - (In) Global row to extract. \param ILength - (In) Length of Values and Indices. \param NumEntries - (Out) Number of nonzero entries extracted. \param Values - (Out) Extracted values for this row. \param Indices - (Out) Extracted global column indices for the corresponding values. \return Integer error code, set to 0 if successful, non-zero if global row is not owned by calling process
or if the number of entries in this row exceed the Length parameter.
Definition at line 1395 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ExtractGlobalRowCopy | ( | long long | GlobalRow, |
int | Length, | ||
int & | NumEntries, | ||
double * | Values, | ||
long long * | Indices | ||
) | const |
Definition at line 1405 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1416 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1448 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ExtractGlobalRowCopy | ( | long long | GlobalRow, |
int | Length, | ||
int & | NumEntries, | ||
double * | Values | ||
) | const |
Definition at line 1457 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1467 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1487 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1558 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ExtractGlobalRowView | ( | long long | GlobalRow, |
int & | NumEntries, | ||
double *& | Values, | ||
long long *& | Indices | ||
) | const |
Definition at line 1567 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1577 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1597 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ExtractGlobalRowView | ( | long long | GlobalRow, |
int & | NumEntries, | ||
double *& | Values | ||
) | const |
Definition at line 1606 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1616 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 3028 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::Multiply1 | ( | bool | TransA, |
const Epetra_Vector & | x, | ||
Epetra_Vector & | y | ||
) | const |
Definition at line 4050 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 3107 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::Multiply1 | ( | bool | TransA, |
const Epetra_MultiVector & | X, | ||
Epetra_MultiVector & | Y | ||
) | const |
Definition at line 4169 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1629 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1671 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1716 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1771 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1930 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1816 of file Epetra_CrsMatrix.cpp.
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. |
Definition at line 1873 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1973 of file Epetra_CrsMatrix.cpp.
|
inlinevirtual |
If FillComplete() has been called, this query returns true, otherwise it returns false.
Implements Epetra_RowMatrix.
Definition at line 1007 of file Epetra_CrsMatrix.h.
|
inline |
If OptimizeStorage() has been called, this query returns true, otherwise it returns false.
Definition at line 1010 of file Epetra_CrsMatrix.h.
|
inline |
If matrix indices has not been transformed to local, this query returns true, otherwise it returns false.
Definition at line 1013 of file Epetra_CrsMatrix.h.
|
inline |
If matrix indices has been transformed to local, this query returns true, otherwise it returns false.
Definition at line 1016 of file Epetra_CrsMatrix.h.
|
inline |
If matrix indices are packed into single array (done in OptimizeStorage()) return true, otherwise false.
Definition at line 1019 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.
Implements Epetra_RowMatrix.
Definition at line 1022 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.
Implements Epetra_RowMatrix.
Definition at line 1025 of file Epetra_CrsMatrix.h.
|
inline |
If matrix has no diagonal entries in global index space, this query returns true, otherwise it returns false.
Definition at line 1028 of file Epetra_CrsMatrix.h.
|
virtual |
Returns the infinity norm of the global matrix.
Implements Epetra_RowMatrix.
Definition at line 2013 of file Epetra_CrsMatrix.cpp.
|
virtual |
Returns the one norm of the global matrix.
Implements Epetra_RowMatrix.
Definition at line 2060 of file Epetra_CrsMatrix.cpp.
double Epetra_CrsMatrix::NormFrobenius | ( | ) | const |
Returns the frobenius norm of the global matrix.
Definition at line 2113 of file Epetra_CrsMatrix.cpp.
|
inlinevirtual |
Returns the number of nonzero entries in the global matrix.
Implements Epetra_RowMatrix.
Definition at line 1066 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Implements Epetra_RowMatrix.
Definition at line 1072 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of global matrix rows.
Implements Epetra_RowMatrix.
Definition at line 1076 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Implements Epetra_RowMatrix.
Definition at line 1082 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of global matrix columns.
Implements Epetra_RowMatrix.
Definition at line 1086 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Implements Epetra_RowMatrix.
Definition at line 1092 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.
Implements Epetra_RowMatrix.
Definition at line 1096 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Implements Epetra_RowMatrix.
Definition at line 1102 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of nonzero entries in the calling processor's portion of the matrix.
Implements Epetra_RowMatrix.
Definition at line 1105 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of matrix rows owned by the calling processor.
Implements Epetra_RowMatrix.
Definition at line 1108 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
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.
Definition at line 1115 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.
Implements Epetra_RowMatrix.
Definition at line 1121 of file Epetra_CrsMatrix.h.
|
inline |
Returns the current number of nonzero entries in specified global row on this processor.
Definition at line 1124 of file Epetra_CrsMatrix.h.
|
inline |
Returns the allocated number of nonzero entries in specified global row on this processor.
Definition at line 1127 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the maximum number of nonzero entries across all rows on this processor.
Implements Epetra_RowMatrix.
Definition at line 1133 of file Epetra_CrsMatrix.h.
|
inline |
Returns the maximum number of nonzero entries across all rows on all processors.
Definition at line 1139 of file Epetra_CrsMatrix.h.
|
inline |
Returns the current number of nonzero entries in specified local row on this processor.
Definition at line 1142 of file Epetra_CrsMatrix.h.
|
inline |
Returns the allocated number of nonzero entries in specified local row on this processor.
Definition at line 1145 of file Epetra_CrsMatrix.h.
|
inline |
Returns the index base for row and column indices for this graph.
Index base for this map.
Definition at line 1150 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1156 of file Epetra_CrsMatrix.h.
|
inline |
Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable.
Definition at line 1160 of file Epetra_CrsMatrix.h.
|
inline |
Returns a reference to the Epetra_CrsGraph object associated with this matrix.
Definition at line 1163 of file Epetra_CrsMatrix.h.
|
inline |
Returns the Epetra_Map object associated with the rows of this matrix.
Definition at line 1166 of file Epetra_CrsMatrix.h.
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.
Definition at line 437 of file Epetra_CrsMatrix.cpp.
|
inline |
Returns true if we have a well-defined ColMap, and returns false otherwise.
Definition at line 1182 of file Epetra_CrsMatrix.h.
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 no entries have been inserted into the matrix (both IndicesAreLocal() and IndicesAreGlobal() are false) or currentmap->PointSameAs(newmap) is true. This is a collective function. Returns 0 if map is replaced, -1 if not.
Definition at line 454 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ReplaceDomainMapAndImporter | ( | const Epetra_Map & | NewDomainMap, |
const Epetra_Import * | NewImporter | ||
) |
Replaces the current DomainMap & Importer with the user-specified map object.
Replaces the current DomainMap and Importer with the user-specified map object, but only if the matrix has been FillCompleted, Importer's TargetMap matches the ColMap and Importer's SourceMap matches the DomainMap (assuming the importer isn't null). If an Importer is passed in, Epetra_CrsMatrix will copy it.
Returns 0 if map/importer is replaced, -1 if not.
Definition at line 471 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::RemoveEmptyProcessesInPlace | ( | const Epetra_BlockMap * | NewMap | ) |
Remove processes owning zero rows from the Maps and their communicator.
Remove processes owning zero rows from the Maps and their communicator.
newMap | [in] This must be the result of calling the removeEmptyProcesses() method on the row Map. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes. |
Definition at line 476 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1230 of file Epetra_CrsMatrix.h.
|
inline |
Returns the Epetra_Map object associated with the domain of this matrix operator.
Definition at line 1236 of file Epetra_CrsMatrix.h.
|
inline |
Returns the Epetra_Map object associated with the range of this matrix operator.
Definition at line 1242 of file Epetra_CrsMatrix.h.
|
inline |
Returns the Epetra_Import object that contains the import operations for distributed operations.
Definition at line 1245 of file Epetra_CrsMatrix.h.
|
inline |
Returns the Epetra_Export object that contains the export operations for distributed operations.
Definition at line 1248 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns a pointer to the Epetra_Comm communicator associated with this matrix.
Implements Epetra_Operator.
Definition at line 1251 of file Epetra_CrsMatrix.h.
|
inline |
Returns the local row index for given global row index, returns -1 if no local row for this global row.
Definition at line 1258 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1262 of file Epetra_CrsMatrix.h.
|
inline |
Returns the global row index for give local row index, returns IndexBase-1 if we don't have this local row.
Definition at line 1273 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1279 of file Epetra_CrsMatrix.h.
|
inline |
Returns the local column index for given global column index, returns -1 if no local column for this global column.
Definition at line 1286 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1289 of file Epetra_CrsMatrix.h.
|
inline |
Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column.
Definition at line 1297 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1303 of file Epetra_CrsMatrix.h.
|
inline |
Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 1307 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1310 of file Epetra_CrsMatrix.h.
|
inline |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 1314 of file Epetra_CrsMatrix.h.
|
inline |
Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 1321 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1324 of file Epetra_CrsMatrix.h.
|
inline |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 1331 of file Epetra_CrsMatrix.h.
|
inline |
Returns true of GID is owned by the calling processor, otherwise it returns false.
Definition at line 1335 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1338 of file Epetra_CrsMatrix.h.
|
virtual |
Print method.
Reimplemented from Epetra_DistObject.
Definition at line 2891 of file Epetra_CrsMatrix.cpp.
|
inlinevirtual |
Returns a character string describing the operator.
Reimplemented from Epetra_Object.
Definition at line 1354 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
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.
Definition at line 1365 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
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.
Definition at line 1376 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
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.
Definition at line 1393 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns true because this class can compute an Inf-norm.
Implements Epetra_Operator.
Definition at line 1397 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the current UseTranspose setting.
Implements Epetra_Operator.
Definition at line 1400 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with the domain of this matrix operator.
Implements Epetra_Operator.
Definition at line 1403 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with the range of this matrix operator.
Implements Epetra_Operator.
Definition at line 1410 of file Epetra_CrsMatrix.h.
|
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.
Definition at line 1428 of file Epetra_CrsMatrix.cpp.
|
inlinevirtual |
Map() method inherited from Epetra_DistObject.
Implements Epetra_SrcDistObject.
Definition at line 1433 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with the rows of this matrix.
Implements Epetra_RowMatrix.
Definition at line 1436 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the Epetra_Map object associated with columns of this matrix.
Implements Epetra_RowMatrix.
Definition at line 1439 of file Epetra_CrsMatrix.h.
|
inlinevirtual |
Returns the Epetra_Import object that contains the import operations for distributed operations.
Implements Epetra_RowMatrix.
Definition at line 1442 of file Epetra_CrsMatrix.h.
|
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. |
Definition at line 1455 of file Epetra_CrsMatrix.h.
|
inline |
Definition at line 1458 of file Epetra_CrsMatrix.h.
|
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. |
Definition at line 1479 of file Epetra_CrsMatrix.h.
Epetra_IntSerialDenseVector & Epetra_CrsMatrix::ExpertExtractIndexOffset | ( | ) |
Returns a reference to the Epetra_IntSerialDenseVector used to hold the local IndexOffsets (CRS rowptr)
Definition at line 4588 of file Epetra_CrsMatrix.cpp.
Epetra_IntSerialDenseVector & Epetra_CrsMatrix::ExpertExtractIndices | ( | ) |
Returns a reference to the Epetra_IntSerialDenseVector used to hold the local All_Indices (CRS colind)
Definition at line 4593 of file Epetra_CrsMatrix.cpp.
|
inline |
Returns a reference to the double* used to hold the values array.
Definition at line 1505 of file Epetra_CrsMatrix.h.
int Epetra_CrsMatrix::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.
Performs a lightweight FillComplete on an object that already has filled IndexOffsets, All_Indices and All_Values. This routine is needed to support the EpetraExt::MatrixMatrix::Multiply and should not be called by users.
Definition at line 4607 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::ExpertMakeUniqueCrsGraphData | ( | ) |
Makes sure this matrix has a unique CrsGraphData object.
This routine is needed to support the EpetraExt::MatrixMatrix::Multiply and should not be called by users.
Definition at line 4598 of file Epetra_CrsMatrix.cpp.
|
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.
Definition at line 1534 of file Epetra_CrsMatrix.h.
|
inline |
Use ColMap() instead.
Definition at line 1541 of file Epetra_CrsMatrix.h.
int Epetra_CrsMatrix::TransformToLocal | ( | ) |
Use FillComplete() instead.
Definition at line 1187 of file Epetra_CrsMatrix.cpp.
int Epetra_CrsMatrix::TransformToLocal | ( | const Epetra_Map * | DomainMap, |
const Epetra_Map * | RangeMap | ||
) |
Use FillComplete(const Epetra_Map& DomainMap, const Epetra_Map& RangeMap) instead.
Definition at line 1193 of file Epetra_CrsMatrix.cpp.
|
inlineprotected |
Definition at line 1553 of file Epetra_CrsMatrix.h.
|
inlineprotected |
Definition at line 1554 of file Epetra_CrsMatrix.h.
|
inlineprotected |
Definition at line 1555 of file Epetra_CrsMatrix.h.
|
inlineprotected |
Definition at line 1558 of file Epetra_CrsMatrix.h.
|
inlineprotected |
Definition at line 1561 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 319 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 335 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 791 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 665 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 802 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 676 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 814 of file Epetra_CrsMatrix.cpp.
|
protected |
|
protected |
Definition at line 900 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 1108 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3287 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3301 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3315 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3415 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3486 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3587 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3676 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 3807 of file Epetra_CrsMatrix.cpp.
|
inlineprotected |
Definition at line 1593 of file Epetra_CrsMatrix.h.
|
protectedvirtual |
Allows the source and target (this) objects to be compared for compatibility, return nonzero if not.
Implements Epetra_DistObject.
Definition at line 2151 of file Epetra_CrsMatrix.cpp.
|
protectedvirtual |
Perform ID copies and permutations that are on processor.
Implements Epetra_DistObject.
Definition at line 2163 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 2417 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 2596 of file Epetra_CrsMatrix.cpp.
|
protectedvirtual |
Perform any packing or preparation required for call to DoTransfer().
Implements Epetra_DistObject.
Definition at line 2625 of file Epetra_CrsMatrix.cpp.
|
protectedvirtual |
Perform any unpacking and combining after call to DoTransfer().
Implements Epetra_DistObject.
Definition at line 2857 of file Epetra_CrsMatrix.cpp.
|
protected |
Sort column entries, row-by-row, in ascending order.
Definition at line 1199 of file Epetra_CrsMatrix.cpp.
|
inlineprotected |
If SortEntries() has been called, this query returns true, otherwise it returns false.
Definition at line 1643 of file Epetra_CrsMatrix.h.
|
protected |
Add entries that have the same column index. Remove redundant entries from list.
Definition at line 1241 of file Epetra_CrsMatrix.cpp.
|
inlineprotected |
If MergeRedundantEntries() has been called, this query returns true, otherwise it returns false.
Definition at line 1649 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 394 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 4305 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 4439 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 523 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 564 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 647 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 688 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 823 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 933 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 1382 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 1439 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 1547 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 1588 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 2195 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 2447 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 2764 of file Epetra_CrsMatrix.cpp.
|
private |
Definition at line 4803 of file Epetra_CrsMatrix.cpp.
void Epetra_CrsMatrix::FusedImport | ( | const Epetra_CrsMatrix & | SourceMatrix, |
const Epetra_Import & | RowImporter, | ||
const Epetra_Map * | DomainMap, | ||
const Epetra_Map * | RangeMap, | ||
bool | RestrictCommunicator | ||
) |
Definition at line 5197 of file Epetra_CrsMatrix.cpp.
void Epetra_CrsMatrix::FusedExport | ( | const Epetra_CrsMatrix & | SourceMatrix, |
const Epetra_Export & | RowExporter, | ||
const Epetra_Map * | DomainMap, | ||
const Epetra_Map * | RangeMap, | ||
bool | RestrictCommunicator | ||
) |
Definition at line 5206 of file Epetra_CrsMatrix.cpp.
void Epetra_CrsMatrix::FusedImport | ( | const Epetra_CrsMatrix & | SourceMatrix, |
const Epetra_Import & | RowImporter, | ||
const Epetra_Import * | DomainImporter, | ||
const Epetra_Map * | DomainMap, | ||
const Epetra_Map * | RangeMap, | ||
bool | RestrictCommunicator | ||
) |
Definition at line 5214 of file Epetra_CrsMatrix.cpp.
void Epetra_CrsMatrix::FusedExport | ( | const Epetra_CrsMatrix & | SourceMatrix, |
const Epetra_Export & | RowExporter, | ||
const Epetra_Export * | DomainExporter, | ||
const Epetra_Map * | DomainMap, | ||
const Epetra_Map * | RangeMap, | ||
bool | RestrictCommunicator | ||
) |
Definition at line 5224 of file Epetra_CrsMatrix.cpp.
|
protected |
Definition at line 1653 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1654 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1655 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1656 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1657 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1658 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1659 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1661 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1662 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1663 of file Epetra_CrsMatrix.h.
|
mutableprotected |
Definition at line 1664 of file Epetra_CrsMatrix.h.
|
mutableprotected |
Definition at line 1665 of file Epetra_CrsMatrix.h.
|
mutableprotected |
Definition at line 1666 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1668 of file Epetra_CrsMatrix.h.
|
mutableprotected |
Definition at line 1669 of file Epetra_CrsMatrix.h.
|
mutableprotected |
Definition at line 1670 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1672 of file Epetra_CrsMatrix.h.
|
protected |
Definition at line 1674 of file Epetra_CrsMatrix.h.