Epetra Package Browser (Single Doxygen Collection) Development
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions
Epetra_CrsMatrix Class Reference

Epetra_CrsMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices. More...

#include <Epetra_CrsMatrix.h>

Inheritance diagram for Epetra_CrsMatrix:
Inheritance graph
[legend]

List of all members.

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)

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.
int CopyAndPermute (const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero)
 Perform ID copies and permutations that are on processor.
int CopyAndPermuteCrsMatrix (const Epetra_CrsMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode)
int CopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode)
int PackAndPrepare (const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor)
 Perform any packing or preparation required for call to DoTransfer().
int UnpackAndCombine (const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor)
 Perform any unpacking and combining after call to DoTransfer().
int SortEntries ()
 Sort column entries, row-by-row, in ascending order.
bool Sorted () const
 If SortEntries() has been called, this query returns true, otherwise it returns false.
int MergeRedundantEntries ()
 Add entries that have the same column index. Remove redundant entries from list.
bool NoRedundancies () const
 If MergeRedundantEntries() has been called, this query returns true, otherwise it returns false.
void DeleteMemory ()

Protected 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_MultiVectorImportVector_
Epetra_MultiVectorExportVector_
Epetra_DataAccess CV_
bool squareFillCompleteCalled_

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 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.
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, int NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constructor with fixed number of indices per row.
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, const int *NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constructor with variable number of indices per row.
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_Map &RowMap, const Epetra_Map &ColMap, int NumEntriesPerRow, bool StaticProfile=false)
 Epetra_CrsMatrix constuctor with fixed number of indices per row.
 Epetra_CrsMatrix (Epetra_DataAccess CV, const Epetra_CrsGraph &Graph)
 Construct a matrix using an existing Epetra_CrsGraph object.
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Import &RowImporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false)
 Epetra CrsMatrix constructor that also fuses Import and FillComplete().
 Epetra_CrsMatrix (const Epetra_CrsMatrix &SourceMatrix, const Epetra_Export &RowExporter, const Epetra_Map *DomainMap=0, const Epetra_Map *RangeMap=0, bool RestrictCommunicator=false)
 Epetra CrsMatrix constructor that also fuses Ex[prt and FillComplete().
 Epetra_CrsMatrix (const Epetra_CrsMatrix &Matrix)
 Copy constructor.
virtual ~Epetra_CrsMatrix ()
 Epetra_CrsMatrix Destructor.

Insertion/Replace/SumInto methods

Epetra_CrsMatrixoperator= (const Epetra_CrsMatrix &src)
 Assignment operator.
int PutScalar (double ScalarConstant)
 Initialize all values in the matrix with constant value.
int Scale (double ScalarConstant)
 Multiply all values in the matrix by a constant value (in place: A <- ScalarConstant * A).
virtual int InsertGlobalValues (int GlobalRow, int NumEntries, const double *Values, const int *Indices)
 Insert a list of elements in a given global row of the matrix.
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.
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.
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.
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.
int InsertMyValues (int MyRow, int NumEntries, double *Values, int *Indices)
 Insert a list of elements in a given local row of the matrix.
int ReplaceMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices)
 Replace current values with this list of entries for a given local row of the matrix.
int SumIntoMyValues (int MyRow, int NumEntries, const double *Values, const int *Indices)
 Add this list of entries to existing values for a given local row of the matrix.
int ReplaceDiagonalValues (const Epetra_Vector &Diagonal)
 Replaces diagonal values of the matrix with those in the user-provided vector.

Transformation methods

int FillComplete (bool OptimizeDataStorage=true)
 Signal that data entry is complete. Perform transformations to local index space.
int FillComplete (const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, bool OptimizeDataStorage=true)
 Signal that data entry is complete. Perform transformations to local index space.
int OptimizeStorage ()
 Make consecutive row index sections contiguous, minimize internal storage used for constructing graph.
int MakeDataContiguous ()
 Eliminates memory that is used for construction. Make consecutive row index sections contiguous.

Extraction methods

int ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified global row in user-provided arrays.
int ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values, long long *Indices) const
int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified local row in user-provided arrays.
int ExtractGlobalRowCopy (int GlobalRow, int Length, int &NumEntries, double *Values) const
 Returns a copy of the specified global row values in user-provided array.
int ExtractGlobalRowCopy (long long GlobalRow, int Length, int &NumEntries, double *Values) const
int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values) const
 Returns a copy of the specified local row values in user-provided array.
int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 Returns a copy of the main diagonal in a user-provided vector.
int ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values, int *&Indices) const
 Returns a view of the specified global row values via pointers to internal data.
int ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values, long long *&Indices) const
int ExtractMyRowView (int MyRow, int &NumEntries, double *&Values, int *&Indices) const
 Returns a view of the specified local row values via pointers to internal data.
int ExtractGlobalRowView (int GlobalRow, int &NumEntries, double *&Values) const
 Returns a view of the specified global row values via pointers to internal data.
int ExtractGlobalRowView (long long GlobalRow, int &NumEntries, double *&Values) const
int ExtractMyRowView (int MyRow, int &NumEntries, double *&Values) const
 Returns a view of the specified local row values via pointers to internal data.

Computational methods

int Multiply (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
 Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_Vector x in y.
int Multiply1 (bool TransA, const Epetra_Vector &x, Epetra_Vector &y) const
int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y.
int Multiply1 (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
int Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_Vector &x, Epetra_Vector &y) const
 Returns the result of a local solve using the Epetra_CrsMatrix on a Epetra_Vector x in y.
int Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a local solve using the Epetra_CrsMatrix a Epetra_MultiVector X in Y.
int InvRowSums (Epetra_Vector &x) const
 Computes the inverse of the sum of absolute values of the rows of the Epetra_CrsMatrix, results returned in x.
int InvRowMaxs (Epetra_Vector &x) const
 Computes the inverse of the max of absolute values of the rows of the Epetra_CrsMatrix, results returned in x.
int LeftScale (const Epetra_Vector &x)
 Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x.
int InvColSums (Epetra_Vector &x) const
 Computes the inverse of the sum of absolute values of the columns of the Epetra_CrsMatrix, results returned in x.
int InvColMaxs (Epetra_Vector &x) const
 Computes the max of absolute values of the columns of the Epetra_CrsMatrix, results returned in x.
int RightScale (const Epetra_Vector &x)
 Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x.

Matrix Properties Query Methods

bool Filled () const
 If FillComplete() has been called, this query returns true, otherwise it returns false.
bool StorageOptimized () const
 If OptimizeStorage() has been called, this query returns true, otherwise it returns false.
bool IndicesAreGlobal () const
 If matrix indices has not been transformed to local, this query returns true, otherwise it returns false.
bool IndicesAreLocal () const
 If matrix indices has been transformed to local, this query returns true, otherwise it returns false.
bool IndicesAreContiguous () const
 If matrix indices are packed into single array (done in OptimizeStorage()) return true, otherwise false.
bool LowerTriangular () const
 If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.
bool UpperTriangular () const
 If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.
bool NoDiagonal () const
 If matrix has no diagonal entries in global index space, this query returns true, otherwise it returns false.

Attribute access functions

double NormInf () const
 Returns the infinity norm of the global matrix.
double NormOne () const
 Returns the one norm of the global matrix.
double NormFrobenius () const
 Returns the frobenius norm of the global matrix.
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global matrix.
long long NumGlobalNonzeros64 () const
int NumGlobalRows () const
 Returns the number of global matrix rows.
long long NumGlobalRows64 () const
int NumGlobalCols () const
 Returns the number of global matrix columns.
long long NumGlobalCols64 () const
int NumGlobalDiagonals () const
 Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.
long long NumGlobalDiagonals64 () const
int NumMyNonzeros () const
 Returns the number of nonzero entries in the calling processor's portion of the matrix.
int NumMyRows () const
 Returns the number of matrix rows owned by the calling processor.
int NumMyCols () const
 Returns the number of entries in the set of column-indices that appear on this processor.
int NumMyDiagonals () const
 Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.
int NumGlobalEntries (long long Row) const
 Returns the current number of nonzero entries in specified global row on this processor.
int NumAllocatedGlobalEntries (int Row) const
 Returns the allocated number of nonzero entries in specified global row on this processor.
int MaxNumEntries () const
 Returns the maximum number of nonzero entries across all rows on this processor.
int GlobalMaxNumEntries () const
 Returns the maximum number of nonzero entries across all rows on all processors.
int NumMyEntries (int Row) const
 Returns the current number of nonzero entries in specified local row on this processor.
int NumAllocatedMyEntries (int Row) const
 Returns the allocated number of nonzero entries in specified local row on this processor.
int IndexBase () const
 Returns the index base for row and column indices for this graph.
long long IndexBase64 () const
bool StaticGraph ()
 Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable.
const Epetra_CrsGraphGraph () const
 Returns a reference to the Epetra_CrsGraph object associated with this matrix.
const Epetra_MapRowMap () const
 Returns the Epetra_Map object associated with the rows of this matrix.
int ReplaceRowMap (const Epetra_BlockMap &newmap)
 Replaces the current RowMap with the user-specified map object.
bool HaveColMap () const
 Returns true if we have a well-defined ColMap, and returns false otherwise.
int ReplaceColMap (const Epetra_BlockMap &newmap)
 Replaces the current ColMap with the user-specified map object.
int ReplaceDomainMapAndImporter (const Epetra_Map &NewDomainMap, const Epetra_Import *NewImporter)
 Replaces the current DomainMap & Importer with the user-specified map object.
int RemoveEmptyProcessesInPlace (const Epetra_BlockMap *NewMap)
 Remove processes owning zero rows from the Maps and their communicator.
const Epetra_MapColMap () const
 Returns the Epetra_Map object that describes the set of column-indices that appear in each processor's locally owned matrix rows.
const Epetra_MapDomainMap () const
 Returns the Epetra_Map object associated with the domain of this matrix operator.
const Epetra_MapRangeMap () const
 Returns the Epetra_Map object associated with the range of this matrix operator.
const Epetra_ImportImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations.
const Epetra_ExportExporter () const
 Returns the Epetra_Export object that contains the export operations for distributed operations.
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix.

Local/Global ID methods

int LRID (int GRID_in) const
 Returns the local row index for given global row index, returns -1 if no local row for this global row.
int LRID (long long GRID_in) const
int GRID (int LRID_in) const
 Returns the global row index for give local row index, returns IndexBase-1 if we don't have this local row.
long long GRID64 (int LRID_in) const
int LCID (int GCID_in) const
 Returns the local column index for given global column index, returns -1 if no local column for this global column.
int LCID (long long GCID_in) const
int GCID (int LCID_in) const
 Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column.
long long GCID64 (int LCID_in) const
bool MyGRID (int GRID_in) const
 Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false.
bool MyGRID (long long GRID_in) const
bool MyLRID (int LRID_in) const
 Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
bool MyGCID (int GCID_in) const
 Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false.
bool MyGCID (long long GCID_in) const
bool MyLCID (int LCID_in) const
 Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
bool MyGlobalRow (int GID) const
 Returns true of GID is owned by the calling processor, otherwise it returns false.
bool MyGlobalRow (long long GID) const

I/O Methods

virtual void Print (std::ostream &os) const
 Print method.

Additional methods required to support the Epetra_Operator interface

const char * Label () const
 Returns a character string describing the operator.
int SetUseTranspose (bool UseTranspose_in)
 If set true, transpose of this operator will be applied.
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.
bool HasNormInf () const
 Returns true because this class can compute an Inf-norm.
bool UseTranspose () const
 Returns the current UseTranspose setting.
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this matrix operator.
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this matrix operator.

Additional methods required to implement Epetra_RowMatrix interface

int NumMyRowEntries (int MyRow, int &NumEntries) const
 Return the current number of values stored for the specified local row.
const Epetra_BlockMapMap () const
 Map() method inherited from Epetra_DistObject.
const Epetra_MapRowMatrixRowMap () const
 Returns the Epetra_Map object associated with the rows of this matrix.
const Epetra_MapRowMatrixColMap () const
 Returns the Epetra_Map object associated with columns of this matrix.
const Epetra_ImportRowMatrixImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations.

Inlined Operator Methods

double * operator[] (int Loc)
 Inlined bracket operator for fast access to data. (Const and Non-const versions)
double * operator[] (int Loc) const

Expert-only methods: These methods are intended for experts only and have some risk of changing in the future, since they rely on underlying data structure assumptions

int ExtractCrsDataPointers (int *&IndexOffset, int *&Indices, double *&Values_in) const
 Returns internal data pointers associated with Crs matrix format.
Epetra_IntSerialDenseVectorExpertExtractIndexOffset ()
 Returns a reference to the Epetra_IntSerialDenseVector used to hold the local IndexOffsets (CRS rowptr)
Epetra_IntSerialDenseVectorExpertExtractIndices ()
 Returns a reference to the Epetra_IntSerialDenseVector used to hold the local All_Indices (CRS colind)
double *& ExpertExtractValues ()
 Returns a reference to the double* used to hold the values array.
int ExpertStaticFillComplete (const Epetra_Map &DomainMap, const Epetra_Map &RangeMap, const Epetra_Import *Importer=0, const Epetra_Export *Exporter=0, int NumMyDiagonals=-1)
 Performs a FillComplete on an object that aready has filled CRS data.
int ExpertMakeUniqueCrsGraphData ()
 Makes sure this matrix has a unique CrsGraphData object.
int SortGhostsAssociatedWithEachProcessor (bool Flag)
 Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order.

Deprecated methods: These methods still work, but will be removed in a future version

const Epetra_MapImportMap () const
 Use ColMap() instead.
int TransformToLocal ()
 Use FillComplete() instead.
int TransformToLocal (const Epetra_Map *DomainMap, const Epetra_Map *RangeMap)
 Use FillComplete(const Epetra_Map& DomainMap, const Epetra_Map& RangeMap) instead.

Detailed Description

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:

  1. Create Epetra_CrsMatrix instance, including storage, via one of the constructors:
    • Constructor that accepts one Epetra_Map object, a row-map defining the distribution of matrix rows.
    • Constructor that accepts two Epetra_Map objects. (The second map is a column-map, and describes the set of column-indices that appear in each processor's portion of the matrix. Generally these are overlapping sets -- column-indices may appear on more than one processor.)
    • Constructor that accepts an Epetra_CrsGraph object, defining the non-zero structure of the matrix.
    Note that the constructors which accept Epetra_Map arguments also accept an argument that gives an estimate of the number of nonzeros per row. This allows storage to be pre-allocated and can improve the performance of the data input methods. The estimate need not be accurate, as additional storage is allocated automatically when needed. However, a more accurate estimate helps performance by reducing the amount of extra memory allocation.
  2. Enter values via one or more Insert/Replace/SumInto functions.
  3. Complete construction by calling FillComplete.

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.

Warning:
A Epetra_Map is required for the Epetra_CrsMatrix constructor.

Definition at line 173 of file Epetra_CrsMatrix.h.


Constructor & Destructor Documentation

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.

Parameters:
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.

Parameters:
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.

Parameters:
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.

Parameters:
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.

Parameters:
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

Parameters:
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 5049 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

Parameters:
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 5067 of file Epetra_CrsMatrix.cpp.

Epetra_CrsMatrix::Epetra_CrsMatrix ( const Epetra_CrsMatrix Matrix)

Copy constructor.

Definition at line 234 of file Epetra_CrsMatrix.cpp.

Epetra_CrsMatrix::~Epetra_CrsMatrix ( ) [virtual]

Epetra_CrsMatrix Destructor.

Definition at line 388 of file Epetra_CrsMatrix.cpp.


Member Function Documentation

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.

Parameters:
ScalarConstant- (In) Value to use.
Returns:
Integer error code, set to 0 if successful.
Precondition:
None.
Postcondition:
All values in this set to ScalarConstant.

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).

Parameters:
ScalarConstant- (In) Value to use.
Returns:
Integer error code, set to 0 if successful.
Precondition:
None.
Postcondition:
All values of this have been multiplied by ScalarConstant.

Definition at line 504 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertGlobalValues ( int  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [virtual]

Insert a list of elements in a given global row of the matrix.

This method is used to construct a matrix for the first time. It cannot be used if the matrix structure has already been fixed (via a call to FillComplete()). If multiple values are inserted for the same matrix entry, the values are initially stored separately, so memory use will grow as a result. However, when FillComplete is called the values will be summed together and the additional memory will be released.

For example, if the values 2.0, 3.0 and 4.0 are all inserted in Row 1, Column 2, extra storage is used to store each of the three values separately. In this way, the insert process does not require any searching and can be faster. However, when FillComplete() is called, the values will be summed together to equal 9.0 and only a single entry will remain in the matrix for Row 1, Column 2.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code will be returned.
Warning:
This method may not be called once FillComplete() has been called.
Precondition:
IndicesAreLocal()==false && IndicesAreContiguous()==false

Reimplemented in Epetra_FECrsMatrix.

Definition at line 540 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertGlobalValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const long long *  Indices 
) [virtual]

Reimplemented in Epetra_FECrsMatrix.

Definition at line 551 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertGlobalValues ( int  GlobalRow,
int  NumEntries,
double *  Values,
int *  Indices 
) [virtual]

Insert a list of elements in a given global row of the matrix.

This method is used to construct a matrix for the first time. It cannot be used if the matrix structure has already been fixed (via a call to FillComplete()). If multiple values are inserted for the same matrix entry, the values are initially stored separately, so memory use will grow as a result. However, when FillComplete is called the values will be summed together and the additional memory will be released.

For example, if the values 2.0, 3.0 and 4.0 are all inserted in Row 1, Column 2, extra storage is used to store each of the three values separately. In this way, the insert process does not require any searching and can be faster. However, when FillComplete() is called, the values will be summed together to equal 9.0 and only a single entry will remain in the matrix for Row 1, Column 2.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code will be returned.
Warning:
This method may not be called once FillComplete() has been called.
Precondition:
IndicesAreLocal()==false && IndicesAreContiguous()==false

Reimplemented in Epetra_FECrsMatrix.

Definition at line 581 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertGlobalValues ( long long  GlobalRow,
int  NumEntries,
double *  Values,
long long *  Indices 
) [virtual]

Reimplemented in Epetra_FECrsMatrix.

Definition at line 592 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::ReplaceGlobalValues ( int  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [virtual]

Replace specified existing values with this list of entries for a given global row of the matrix.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Note that if a value is not already present for the specified location in the matrix, the input value will be ignored and a positive warning code will be returned.
Precondition:
IndicesAreLocal()==false && IndicesAreContiguous()==false

Reimplemented in Epetra_FECrsMatrix.

Definition at line 852 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::ReplaceGlobalValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const long long *  Indices 
) [virtual]

Reimplemented in Epetra_FECrsMatrix.

Definition at line 860 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::SumIntoGlobalValues ( int  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [virtual]

Add this list of entries to existing values for a given global row of the matrix.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Note that if a value is not already present for the specified location in the matrix, the input value will be ignored and a positive warning code will be returned.
Precondition:
IndicesAreLocal()==false && IndicesAreContiguous()==false

Reimplemented in Epetra_FECrsMatrix.

Definition at line 1028 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::SumIntoGlobalValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const long long *  Indices 
) [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.

Parameters:
MyRow- (In) Row number (in local coordinates) to put elements.
NumEntries- (In) Number of entries.
Values- (In) Values to enter.
Indices- (In) Local column indices corresponding to values.
Returns:
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code will be returned.
Precondition:
IndicesAreGlobal()==false && (IndicesAreContiguous()==false || CV_==View)
Postcondition:
The given local row of the matrix has been updated as described above.

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.

Parameters:
MyRow- (In) Row number (in local coordinates) to put elements.
NumEntries- (In) Number of entries.
Values- (In) Values to enter.
Indices- (In) Local column indices corresponding to values.
Returns:
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code will be returned.
Precondition:
IndicesAreGlobal()==false && (IndicesAreContiguous()==false || CV_==View)
Postcondition:
The given local row of the matrix has been updated as described above.

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.

Parameters:
MyRow- (In) Row number (in local coordinates) to put elements.
NumEntries- (In) Number of entries.
Values- (In) Values to enter.
Indices- (In) Local column indices corresponding to values.
Returns:
Integer error code, set to 0 if successful. Note that if a value is not already present for the specified location in the matrix, the input value will be ignored and a positive warning code will be returned.
Precondition:
IndicesAreLocal()==true
Postcondition:
MyRow contains the given list of Values at the given Indices.

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.

Parameters:
MyRow- (In) Row number (in local coordinates) to put elements.
NumEntries- (In) Number of entries.
Values- (In) Values to enter.
Indices- (In) Local column indices corresponding to values.
Returns:
Integer error code, set to 0 if successful. Note that if the allocated length of the row has to be expanded, a positive warning code will be returned.
Precondition:
IndicesAreLocal()==true
Postcondition:
The given Values at the given Indices have been summed into the entries of MyRow.

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.

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

Note:
A call to this method will also call the OptimizeStorage method for the associated Epetra_CrsGraph object. If the storage for this graph has already been optimized this additional call will have no effect.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true.
If CV=View when the graph was constructed, then this method will be effective if the indices of the graph were already contiguous. In this case, the indices are left untouched and internal storage for the graph is minimized.
Postcondition:
StorageOptimized()==true, if successful.
Graph().StorageOptimized()==true, if successful.

Definition at line 1279 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::MakeDataContiguous ( ) [inline]

Eliminates memory that is used for construction. Make consecutive row index sections contiguous.

Definition at line 632 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.

Parameters:
GlobalRow- (In) Global row to extract.
ILength- (In) Length of Values and Indices.
NumEntries- (Out) Number of nonzero entries extracted.
Values- (Out) Extracted values for this row.
Indices- (Out) Extracted global column indices for the corresponding values.
Returns:
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.

int Epetra_CrsMatrix::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const [virtual]

Returns a copy of the specified local row in user-provided arrays.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.
Precondition:
IndicesAreLocal()==true

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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.

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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.

Definition at line 1467 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::ExtractDiagonalCopy ( Epetra_Vector Diagonal) const [virtual]

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

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

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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Returns -1 of row not on this processor. Returns -2 if matrix is not in global form (if FillComplete() has already been called).
Precondition:
IndicesAreGlobal()==true

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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Returns -1 of row not on this processor. Returns -2 if matrix is not in local form (if FillComplete() has not been called).
Precondition:
IndicesAreLocal()==true

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.

Parameters:
GlobalRow- (In) Global row to extract.
NumEntries- (Out) Number of nonzero entries extracted.
Values- (Out) Extracted values for this row.
Returns:
Integer error code, set to 0 if successful.

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.

Parameters:
MyRow- (In) Local row to extract.
NumEntries- (Out) Number of nonzero entries extracted.
Values- (Out) Extracted values for this row.
Returns:
Integer error code, set to 0 if successful.

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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Reimplemented in Epetra_OskiMatrix.

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 4045 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const [virtual]

Returns the result of a Epetra_CrsMatrix multiplied by a Epetra_MultiVector X in Y.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

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 4164 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.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Reimplemented in Epetra_OskiMatrix.

Definition at line 1629 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::Solve ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const [virtual]

Returns the result of a local solve using the Epetra_CrsMatrix a Epetra_MultiVector X in Y.

This method solves a triangular system of equations asynchronously on each processor.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Implements Epetra_RowMatrix.

Reimplemented in Epetra_OskiMatrix.

Definition at line 1671 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InvRowSums ( Epetra_Vector x) const [virtual]

Computes the inverse of the sum of absolute values of the rows of the Epetra_CrsMatrix, results returned in x.

The vector x will return such that x[i] will contain the inverse of the sum of the absolute values of the entries in the ith row of the this matrix. Using the resulting vector from this function as input to LeftScale() will make the infinity norm of the resulting matrix exactly 1.

Warning:
The NormInf() method will not properly calculate the infinity norm for a matrix that has entries that are replicated on multiple processors. In this case, if the rows are fully replicated, NormInf() will return a value equal to the maximum number of processors that any individual row of the matrix is replicated on.
Parameters:
x- (Out) An Epetra_Vector containing the inverse of the row sums of the this matrix.
Warning:
When rows are fully replicated on multiple processors, it is assumed that the distribution of x is the same as the rows (RowMap())of this. When multiple processors contain partial sums for individual entries, the distribution of x is assumed to be the same as the RangeMap() of this. When each row of this is uniquely owned, the distribution of x can be that of the RowMap() or the RangeMap().
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

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.

Warning:
This method will not work when multiple processors contain partial sums for individual entries.
Parameters:
x- (Out) An Epetra_Vector containing the inverse of the row maxs of the this matrix.
Warning:
When rows are fully replicated on multiple processors, it is assumed that the distribution of x is the same as the rows (RowMap())of this. When each row of this is uniquely owned, the distribution of x can be that of the RowMap() or the RangeMap().
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Definition at line 1771 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::LeftScale ( const Epetra_Vector x) [virtual]

Scales the Epetra_CrsMatrix on the left with a Epetra_Vector x.

The this matrix will be scaled such that A(i,j) = x(i)*A(i,j) where i denotes the row number of A and j denotes the column number of A.

Parameters:
x- (In) An Epetra_Vector to scale with.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
The matrix will be scaled as described above.

Implements Epetra_RowMatrix.

Definition at line 1930 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InvColSums ( Epetra_Vector x) const [virtual]

Computes the inverse of the sum of absolute values of the columns of the Epetra_CrsMatrix, results returned in x.

The vector x will return such that x[j] will contain the inverse of the sum of the absolute values of the entries in the jth column of the this matrix. Using the resulting vector from this function as input to RightScale() will make the one norm of the resulting matrix exactly 1.

Warning:
The NormOne() method will not properly calculate the one norm for a matrix that has entries that are replicated on multiple processors. In this case, if the columns are fully replicated, NormOne() will return a value equal to the maximum number of processors that any individual column of the matrix is repliated on.
Parameters:
x- (Out) An Epetra_Vector containing the column sums of the this matrix.
Warning:
When columns are fully replicated on multiple processors, it is assumed that the distribution of x is the same as the columns (ColMap()) of this. When multiple processors contain partial sums for entries, the distribution of x is assumed to be the same as the DomainMap() of this. When each column of this is uniquely owned, the distribution of x can be that of the ColMap() or the DomainMap().
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

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.

Warning:
This method will not work when multiple processors contain partial sums for individual entries.
Parameters:
x- (Out) An Epetra_Vector containing the column maxs of the this matrix.
Warning:
When columns are fully replicated on multiple processors, it is assumed that the distribution of x is the same as the columns (ColMap()) of this. When each column of this is uniquely owned, the distribution of x can be that of the ColMap() or the DomainMap().
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Definition at line 1873 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::RightScale ( const Epetra_Vector x) [virtual]

Scales the Epetra_CrsMatrix on the right with a Epetra_Vector x.

The this matrix will be scaled such that A(i,j) = x(j)*A(i,j) where i denotes the global row number of A and j denotes the global column number of A.

Parameters:
x- (In) The Epetra_Vector used for scaling this.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
The matrix will be scaled as described above.

Implements Epetra_RowMatrix.

Definition at line 1973 of file Epetra_CrsMatrix.cpp.

bool Epetra_CrsMatrix::Filled ( ) const [inline, virtual]

If FillComplete() has been called, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 921 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::StorageOptimized ( ) const [inline]

If OptimizeStorage() has been called, this query returns true, otherwise it returns false.

Definition at line 924 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::IndicesAreGlobal ( ) const [inline]

If matrix indices has not been transformed to local, this query returns true, otherwise it returns false.

Definition at line 927 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::IndicesAreLocal ( ) const [inline]

If matrix indices has been transformed to local, this query returns true, otherwise it returns false.

Definition at line 930 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::IndicesAreContiguous ( ) const [inline]

If matrix indices are packed into single array (done in OptimizeStorage()) return true, otherwise false.

Definition at line 933 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::LowerTriangular ( ) const [inline, virtual]

If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 936 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::UpperTriangular ( ) const [inline, virtual]

If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 939 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::NoDiagonal ( ) const [inline]

If matrix has no diagonal entries in global index space, this query returns true, otherwise it returns false.

Definition at line 942 of file Epetra_CrsMatrix.h.

double Epetra_CrsMatrix::NormInf ( ) const [virtual]

Returns the infinity norm of the global matrix.

Implements Epetra_RowMatrix.

Definition at line 2013 of file Epetra_CrsMatrix.cpp.

double Epetra_CrsMatrix::NormOne ( ) const [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.

int Epetra_CrsMatrix::NumGlobalNonzeros ( ) const [inline, virtual]

Returns the number of nonzero entries in the global matrix.

Implements Epetra_RowMatrix.

Definition at line 980 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::NumGlobalNonzeros64 ( ) const [inline, virtual]

Implements Epetra_RowMatrix.

Definition at line 986 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumGlobalRows ( ) const [inline, virtual]

Returns the number of global matrix rows.

Implements Epetra_RowMatrix.

Definition at line 990 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::NumGlobalRows64 ( ) const [inline, virtual]

Implements Epetra_RowMatrix.

Definition at line 996 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumGlobalCols ( ) const [inline, virtual]

Returns the number of global matrix columns.

Implements Epetra_RowMatrix.

Definition at line 1000 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::NumGlobalCols64 ( ) const [inline, virtual]

Implements Epetra_RowMatrix.

Definition at line 1006 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumGlobalDiagonals ( ) const [inline, virtual]

Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.

Implements Epetra_RowMatrix.

Definition at line 1010 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::NumGlobalDiagonals64 ( ) const [inline, virtual]

Implements Epetra_RowMatrix.

Definition at line 1016 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyNonzeros ( ) const [inline, virtual]

Returns the number of nonzero entries in the calling processor's portion of the matrix.

Implements Epetra_RowMatrix.

Definition at line 1019 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyRows ( ) const [inline, virtual]

Returns the number of matrix rows owned by the calling processor.

Implements Epetra_RowMatrix.

Definition at line 1022 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyCols ( ) const [inline, virtual]

Returns the number of entries in the set of column-indices that appear on this processor.

The set of column-indices that appear on this processor is the union of column-indices that appear in all local rows. The size of this set isn't available until FillComplete() has been called.

Precondition:
Filled()==true

Implements Epetra_RowMatrix.

Definition at line 1029 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyDiagonals ( ) const [inline, virtual]

Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.

Precondition:
Filled()==true

Implements Epetra_RowMatrix.

Definition at line 1035 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumGlobalEntries ( long long  Row) const [inline]

Returns the current number of nonzero entries in specified global row on this processor.

Definition at line 1038 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumAllocatedGlobalEntries ( int  Row) const [inline]

Returns the allocated number of nonzero entries in specified global row on this processor.

Definition at line 1041 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::MaxNumEntries ( ) const [inline, virtual]

Returns the maximum number of nonzero entries across all rows on this processor.

Precondition:
Filled()==true

Implements Epetra_RowMatrix.

Definition at line 1047 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::GlobalMaxNumEntries ( ) const [inline]

Returns the maximum number of nonzero entries across all rows on all processors.

Precondition:
Filled()==true

Definition at line 1053 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyEntries ( int  Row) const [inline]

Returns the current number of nonzero entries in specified local row on this processor.

Definition at line 1056 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumAllocatedMyEntries ( int  Row) const [inline]

Returns the allocated number of nonzero entries in specified local row on this processor.

Definition at line 1059 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::IndexBase ( ) const [inline]

Returns the index base for row and column indices for this graph.

Index base for this map.

Definition at line 1064 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::IndexBase64 ( ) const [inline]

Definition at line 1070 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::StaticGraph ( ) [inline]

Returns true if the graph associated with this matrix was pre-constructed and therefore not changeable.

Definition at line 1074 of file Epetra_CrsMatrix.h.

const Epetra_CrsGraph& Epetra_CrsMatrix::Graph ( ) const [inline]

Returns a reference to the Epetra_CrsGraph object associated with this matrix.

Definition at line 1077 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::RowMap ( ) const [inline]

Returns the Epetra_Map object associated with the rows of this matrix.

Definition at line 1080 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.

Precondition:
RowMap().PointSameAs(newmap)==true

Definition at line 437 of file Epetra_CrsMatrix.cpp.

bool Epetra_CrsMatrix::HaveColMap ( ) const [inline]

Returns true if we have a well-defined ColMap, and returns false otherwise.

Precondition:
We have a well-defined ColMap if a) a ColMap was passed in at construction, or b) the MakeColMap function has been called. (Calling either of the FillComplete functions will result in MakeColMap being called.)

Definition at line 1096 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.

Precondition:
(IndicesAreLocal()==false && IndicesAreGlobal()==false) || ColMap().PointSameAs(newmap)==true

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.

Precondition:
(!NewImporter && ColMap().PointSameAs(NewDomainMap)) || (NewImporter && ColMap().PointSameAs(NewImporter->TargetMap()) && NewDomainMap.PointSameAs(NewImporter->SourceMap()))

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.

Warning:
This method is ONLY for use by experts.
We make NO promises of backwards compatibility. This method may change or disappear at any time.
Parameters:
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.

const Epetra_Map& Epetra_CrsMatrix::ColMap ( ) const [inline]

Returns the Epetra_Map object that describes the set of column-indices that appear in each processor's locally owned matrix rows.

Note that if the matrix was constructed with only a row-map, then until FillComplete() is called, this method returns a column-map that is a copy of the row-map. That 'initial' column-map is replaced with a computed column-map (that contains the set of column-indices appearing in each processor's local portion of the matrix) when FillComplete() is called.

Precondition:
HaveColMap()==true

Definition at line 1144 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::DomainMap ( ) const [inline]

Returns the Epetra_Map object associated with the domain of this matrix operator.

Precondition:
Filled()==true

Definition at line 1150 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::RangeMap ( ) const [inline]

Returns the Epetra_Map object associated with the range of this matrix operator.

Precondition:
Filled()==true

Definition at line 1156 of file Epetra_CrsMatrix.h.

const Epetra_Import* Epetra_CrsMatrix::Importer ( ) const [inline]

Returns the Epetra_Import object that contains the import operations for distributed operations.

Definition at line 1159 of file Epetra_CrsMatrix.h.

const Epetra_Export* Epetra_CrsMatrix::Exporter ( ) const [inline]

Returns the Epetra_Export object that contains the export operations for distributed operations.

Definition at line 1162 of file Epetra_CrsMatrix.h.

const Epetra_Comm& Epetra_CrsMatrix::Comm ( ) const [inline]

Returns a pointer to the Epetra_Comm communicator associated with this matrix.

Reimplemented from Epetra_DistObject.

Definition at line 1165 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::LRID ( int  GRID_in) const [inline]

Returns the local row index for given global row index, returns -1 if no local row for this global row.

Definition at line 1172 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::LRID ( long long  GRID_in) const [inline]

Definition at line 1176 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::GRID ( int  LRID_in) const [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 1187 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::GRID64 ( int  LRID_in) const [inline]

Definition at line 1193 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::LCID ( int  GCID_in) const [inline]

Returns the local column index for given global column index, returns -1 if no local column for this global column.

Precondition:
HaveColMap()==true (If HaveColMap()==false, returns -1)

Definition at line 1200 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::LCID ( long long  GCID_in) const [inline]

Definition at line 1203 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::GCID ( int  LCID_in) const [inline]

Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column.

Precondition:
HaveColMap()==true (If HaveColMap()==false, returns -1)

Definition at line 1211 of file Epetra_CrsMatrix.h.

long long Epetra_CrsMatrix::GCID64 ( int  LCID_in) const [inline]

Definition at line 1217 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGRID ( int  GRID_in) const [inline]

Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false.

Definition at line 1221 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGRID ( long long  GRID_in) const [inline]

Definition at line 1224 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyLRID ( int  LRID_in) const [inline]

Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.

Definition at line 1228 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGCID ( int  GCID_in) const [inline]

Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false.

Precondition:
HaveColMap()==true (If HaveColMap()==false, returns -1)

Definition at line 1235 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGCID ( long long  GCID_in) const [inline]

Definition at line 1238 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyLCID ( int  LCID_in) const [inline]

Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.

Precondition:
HaveColMap()==true (If HaveColMap()==false, returns -1)

Definition at line 1245 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGlobalRow ( int  GID) const [inline]

Returns true of GID is owned by the calling processor, otherwise it returns false.

Definition at line 1249 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::MyGlobalRow ( long long  GID) const [inline]

Definition at line 1252 of file Epetra_CrsMatrix.h.

void Epetra_CrsMatrix::Print ( std::ostream &  os) const [virtual]

Print method.

Reimplemented from Epetra_DistObject.

Definition at line 2891 of file Epetra_CrsMatrix.cpp.

const char* Epetra_CrsMatrix::Label ( ) const [inline, virtual]

Returns a character string describing the operator.

Reimplemented from Epetra_Object.

Definition at line 1268 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::SetUseTranspose ( bool  UseTranspose_in) [inline, virtual]

If set true, transpose of this operator will be applied.

This flag allows the transpose of the given operator to be used implicitly. Setting this flag affects only the Apply() and ApplyInverse() methods. If the implementation of this interface does not support transpose use, this method should return a value of -1.

Parameters:
UseTranspose- (In) If true, multiply by the transpose of operator, otherwise just use operator.
Returns:
Always returns 0.

Implements Epetra_Operator.

Definition at line 1279 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const [inline, virtual]

Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.

Parameters:
X- (In) An Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y-(Out) An Epetra_MultiVector of dimension NumVectors containing result.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Implements Epetra_Operator.

Definition at line 1290 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const [inline, virtual]

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

In this implementation, we use several existing attributes to determine how virtual method ApplyInverse() should call the concrete method Solve(). We pass in the UpperTriangular(), the Epetra_CrsMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.

Parameters:
X- (In) An Epetra_MultiVector of dimension NumVectors to solve for.
Y- (Out) An Epetra_MultiVector of dimension NumVectors containing result.
Returns:
Integer error code, set to 0 if successful.
Precondition:
Filled()==true
Postcondition:
Unchanged.

Implements Epetra_Operator.

Definition at line 1307 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::HasNormInf ( ) const [inline, virtual]

Returns true because this class can compute an Inf-norm.

Implements Epetra_Operator.

Definition at line 1311 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::UseTranspose ( ) const [inline, virtual]

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 1314 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::OperatorDomainMap ( ) const [inline, virtual]

Returns the Epetra_Map object associated with the domain of this matrix operator.

Implements Epetra_Operator.

Definition at line 1317 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::OperatorRangeMap ( ) const [inline, virtual]

Returns the Epetra_Map object associated with the range of this matrix operator.

Implements Epetra_Operator.

Definition at line 1324 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const [virtual]

Return the current number of values stored for the specified local row.

Similar to NumMyEntries() except NumEntries is returned as an argument and error checking is done on the input value MyRow.

Parameters:
MyRow- (In) Local row.
NumEntries- (Out) Number of nonzero values.
Returns:
Integer error code, set to 0 if successful.
Precondition:
None.
Postcondition:
Unchanged.

Implements Epetra_RowMatrix.

Definition at line 1428 of file Epetra_CrsMatrix.cpp.

const Epetra_BlockMap& Epetra_CrsMatrix::Map ( ) const [inline, virtual]

Map() method inherited from Epetra_DistObject.

Reimplemented from Epetra_DistObject.

Definition at line 1347 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::RowMatrixRowMap ( ) const [inline, virtual]

Returns the Epetra_Map object associated with the rows of this matrix.

Implements Epetra_RowMatrix.

Definition at line 1350 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::RowMatrixColMap ( ) const [inline, virtual]

Returns the Epetra_Map object associated with columns of this matrix.

Implements Epetra_RowMatrix.

Definition at line 1353 of file Epetra_CrsMatrix.h.

const Epetra_Import* Epetra_CrsMatrix::RowMatrixImporter ( ) const [inline, virtual]

Returns the Epetra_Import object that contains the import operations for distributed operations.

Implements Epetra_RowMatrix.

Definition at line 1356 of file Epetra_CrsMatrix.h.

double* Epetra_CrsMatrix::operator[] ( int  Loc) [inline]

Inlined bracket operator for fast access to data. (Const and Non-const versions)

No error checking and dangerous for optimization purposes.

Parameters:
Loc- (In) Local row.
Returns:
reference to pointer to locally indexed Loc row in matrix.

Definition at line 1369 of file Epetra_CrsMatrix.h.

double* Epetra_CrsMatrix::operator[] ( int  Loc) const [inline]

Definition at line 1372 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::ExtractCrsDataPointers ( int *&  IndexOffset,
int *&  Indices,
double *&  Values_in 
) const [inline]

Returns internal data pointers associated with Crs matrix format.

Returns data pointers to facilitate optimized code within external packages.

Parameters:
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.
Returns:
Integer error code, set to 0 if successful. Returns -1 if FillComplete has not been performed or Storage has not been Optimized.
Warning:
This method is intended for expert only, its use may require user code modifications in future versions of Epetra.

Definition at line 1393 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)

Warning:
This method is intended for experts only, its use may require user code modifications in future versions of Epetra.

Definition at line 4583 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)

Warning:
This method is intended for experts only, its use may require user code modifications in future versions of Epetra.

Definition at line 4588 of file Epetra_CrsMatrix.cpp.

double*& Epetra_CrsMatrix::ExpertExtractValues ( ) [inline]

Returns a reference to the double* used to hold the values array.

Warning:
This method is intended for experts only, its use may require user code modifications in future versions of Epetra.

Definition at line 1419 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.

Warning:
Epetra_CrsMatrix will assume ownership of the Importer/Exporter you pass in. You should not deallocate it afterwards.
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 4602 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.

Warning:
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 4593 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::SortGhostsAssociatedWithEachProcessor ( bool  Flag) [inline]

Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order.

To be compliant with AztecOO, FillComplete() already locally orders ghostnodes such that information received from processor k has a lower local numbering than information received from processor j if k is less than j. SortGhostsAssociatedWithEachProcessor(True) further forces FillComplete() to locally number all ghostnodes received from processor k in ascending order. That is, the local numbering of b is less than c if the global numbering of b is less than c and if both b and c are owned by the same processor. This is done to be compliant with some limited block features within ML. In particular, some ML features require that a block structure of the matrix be maintained even within the ghost variables. Always returns 0.

Definition at line 1448 of file Epetra_CrsMatrix.h.

const Epetra_Map& Epetra_CrsMatrix::ImportMap ( ) const [inline]

Use ColMap() instead.

Definition at line 1455 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.

bool Epetra_CrsMatrix::Allocated ( ) const [inline, protected]

Definition at line 1467 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::SetAllocated ( bool  Flag) [inline, protected]

Definition at line 1468 of file Epetra_CrsMatrix.h.

double** Epetra_CrsMatrix::Values ( ) const [inline, protected]

Definition at line 1469 of file Epetra_CrsMatrix.h.

double* Epetra_CrsMatrix::All_Values ( ) const [inline, protected]

Definition at line 1472 of file Epetra_CrsMatrix.h.

double* Epetra_CrsMatrix::Values ( int  LocalRow) const [inline, protected]

Definition at line 1475 of file Epetra_CrsMatrix.h.

void Epetra_CrsMatrix::InitializeDefaults ( ) [protected]

Definition at line 319 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::Allocate ( ) [protected]

Definition at line 335 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertValues ( int  LocalRow,
int  NumEntries,
double *  Values,
int *  Indices 
) [protected]

Definition at line 791 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertValues ( int  LocalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [protected]

Definition at line 665 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertValues ( int  LocalRow,
int  NumEntries,
double *  Values,
long long *  Indices 
) [protected]

Definition at line 802 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertValues ( int  LocalRow,
int  NumEntries,
const double *  Values,
const long long *  Indices 
) [protected]

Definition at line 676 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertOffsetValues ( long long  GlobalRow,
int  NumEntries,
double *  Values,
int *  Indices 
) [protected]

Definition at line 814 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::InsertOffsetValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [protected]
int Epetra_CrsMatrix::ReplaceOffsetValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [protected]

Definition at line 900 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::SumIntoOffsetValues ( long long  GlobalRow,
int  NumEntries,
const double *  Values,
const int *  Indices 
) [protected]

Definition at line 1108 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::UpdateImportVector ( int  NumVectors) const [protected]

Definition at line 3287 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::UpdateExportVector ( int  NumVectors) const [protected]

Definition at line 3301 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralMV ( double *  x,
double *  y 
) const [protected]

Definition at line 3315 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralMTV ( double *  x,
double *  y 
) const [protected]

Definition at line 3410 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralMM ( double **  X,
int  LDX,
double **  Y,
int  LDY,
int  NumVectors 
) const [protected]

Definition at line 3481 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralMTM ( double **  X,
int  LDX,
double **  Y,
int  LDY,
int  NumVectors 
) const [protected]

Definition at line 3582 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralSV ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
double *  x,
double *  y 
) const [protected]

Definition at line 3671 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::GeneralSM ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
double **  X,
int  LDX,
double **  Y,
int  LDY,
int  NumVectors 
) const [protected]

Definition at line 3802 of file Epetra_CrsMatrix.cpp.

void Epetra_CrsMatrix::SetStaticGraph ( bool  Flag) [inline, protected]

Definition at line 1507 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::CheckSizes ( const Epetra_SrcDistObject Source) [protected, virtual]

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.

int Epetra_CrsMatrix::CopyAndPermute ( const Epetra_SrcDistObject Source,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode = Zero 
) [protected, virtual]

Perform ID copies and permutations that are on processor.

Implements Epetra_DistObject.

Definition at line 2163 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::CopyAndPermuteCrsMatrix ( const Epetra_CrsMatrix A,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode 
) [protected]

Definition at line 2417 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::CopyAndPermuteRowMatrix ( const Epetra_RowMatrix A,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode 
) [protected]

Definition at line 2596 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::PackAndPrepare ( const Epetra_SrcDistObject Source,
int  NumExportIDs,
int *  ExportLIDs,
int &  LenExports,
char *&  Exports,
int &  SizeOfPacket,
int *  Sizes,
bool &  VarSizes,
Epetra_Distributor Distor 
) [protected, virtual]

Perform any packing or preparation required for call to DoTransfer().

Implements Epetra_DistObject.

Definition at line 2625 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::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 
) [protected, virtual]

Perform any unpacking and combining after call to DoTransfer().

Implements Epetra_DistObject.

Definition at line 2857 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::SortEntries ( ) [protected]

Sort column entries, row-by-row, in ascending order.

Definition at line 1199 of file Epetra_CrsMatrix.cpp.

bool Epetra_CrsMatrix::Sorted ( ) const [inline, protected]

If SortEntries() has been called, this query returns true, otherwise it returns false.

Definition at line 1557 of file Epetra_CrsMatrix.h.

int Epetra_CrsMatrix::MergeRedundantEntries ( ) [protected]

Add entries that have the same column index. Remove redundant entries from list.

Definition at line 1241 of file Epetra_CrsMatrix.cpp.

bool Epetra_CrsMatrix::NoRedundancies ( ) const [inline, protected]

If MergeRedundantEntries() has been called, this query returns true, otherwise it returns false.

Definition at line 1563 of file Epetra_CrsMatrix.h.

void Epetra_CrsMatrix::DeleteMemory ( ) [protected]

Reimplemented in Epetra_FECrsMatrix.

Definition at line 394 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::Solve1 ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_Vector x,
Epetra_Vector y 
) const [private]

Definition at line 4300 of file Epetra_CrsMatrix.cpp.

int Epetra_CrsMatrix::Solve1 ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const [private]

Definition at line 4434 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TInsertGlobalValues ( int_type  Row,
int  NumEntries,
const double *  values,
const int_type *  Indices 
) [private]

Definition at line 523 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TInsertGlobalValues ( int_type  Row,
int  NumEntries,
double *  values,
int_type *  Indices 
) [private]

Definition at line 564 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::InsertValues ( int  Row,
int  NumEntries,
const double *  values,
const int_type *  Indices 
) [private]

Definition at line 647 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::InsertValues ( int  Row,
int  NumEntries,
double *  values,
int_type *  Indices 
) [private]

Definition at line 688 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TReplaceGlobalValues ( int_type  Row,
int  NumEntries,
const double *  srcValues,
const int_type *  Indices 
) [private]

Definition at line 823 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TSumIntoGlobalValues ( int_type  Row,
int  NumEntries,
const double *  srcValues,
const int_type *  Indices 
) [private]

Definition at line 933 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::ExtractGlobalRowCopy ( int_type  Row,
int  Length,
int &  NumEntries,
double *  values,
int_type *  Indices 
) const [private]

Definition at line 1382 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::ExtractGlobalRowCopy ( int_type  Row,
int  Length,
int &  NumEntries,
double *  values 
) const [private]

Definition at line 1439 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::ExtractGlobalRowView ( int_type  Row,
int &  NumEntries,
double *&  values,
int_type *&  Indices 
) const [private]

Definition at line 1547 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::ExtractGlobalRowView ( int_type  Row,
int &  NumEntries,
double *&  values 
) const [private]

Definition at line 1588 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TCopyAndPermuteCrsMatrix ( const Epetra_CrsMatrix A,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode 
) [private]

Definition at line 2195 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::TCopyAndPermuteRowMatrix ( const Epetra_RowMatrix A,
int  NumSameIDs,
int  NumPermuteIDs,
int *  PermuteToLIDs,
int *  PermuteFromLIDs,
const Epetra_OffsetIndex Indexor,
Epetra_CombineMode  CombineMode 
) [private]

Definition at line 2447 of file Epetra_CrsMatrix.cpp.

template<typename int_type >
int Epetra_CrsMatrix::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 
) [private]

Definition at line 2764 of file Epetra_CrsMatrix.cpp.

template<class TransferType >
void Epetra_CrsMatrix::FusedTransfer ( const Epetra_CrsMatrix SourceMatrix,
const TransferType &  RowTransfer,
const Epetra_Map DomainMap,
const Epetra_Map RangeMap,
bool  RestrictCommunicator 
) [private]

Definition at line 4787 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 5085 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 5094 of file Epetra_CrsMatrix.cpp.


Member Data Documentation

Definition at line 1567 of file Epetra_CrsMatrix.h.

bool Epetra_CrsMatrix::Allocated_ [protected]

Definition at line 1568 of file Epetra_CrsMatrix.h.

Definition at line 1569 of file Epetra_CrsMatrix.h.

Definition at line 1570 of file Epetra_CrsMatrix.h.

Definition at line 1571 of file Epetra_CrsMatrix.h.

Definition at line 1572 of file Epetra_CrsMatrix.h.

Definition at line 1573 of file Epetra_CrsMatrix.h.

double** Epetra_CrsMatrix::Values_ [protected]

Definition at line 1575 of file Epetra_CrsMatrix.h.

Definition at line 1576 of file Epetra_CrsMatrix.h.

double* Epetra_CrsMatrix::All_Values_ [protected]

Definition at line 1577 of file Epetra_CrsMatrix.h.

double Epetra_CrsMatrix::NormInf_ [mutable, protected]

Definition at line 1578 of file Epetra_CrsMatrix.h.

double Epetra_CrsMatrix::NormOne_ [mutable, protected]

Definition at line 1579 of file Epetra_CrsMatrix.h.

double Epetra_CrsMatrix::NormFrob_ [mutable, protected]

Definition at line 1580 of file Epetra_CrsMatrix.h.

Definition at line 1582 of file Epetra_CrsMatrix.h.

Definition at line 1583 of file Epetra_CrsMatrix.h.

Definition at line 1584 of file Epetra_CrsMatrix.h.

Definition at line 1586 of file Epetra_CrsMatrix.h.

Definition at line 1588 of file Epetra_CrsMatrix.h.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines