Epetra Package Browser (Single Doxygen Collection)
Development
|
Epetra_CrsGraph: A class for constructing and using sparse compressed row graphs. More...
#include <Epetra_CrsGraph.h>
Public Member Functions | |
Epetra_CrsGraph & | operator= (const Epetra_CrsGraph &Source) |
Assignment operator. More... | |
Public Member Functions inherited from Epetra_DistObject | |
Epetra_DistObject (const Epetra_BlockMap &Map) | |
Basic Epetra_DistObject constuctor. More... | |
Epetra_DistObject (const Epetra_BlockMap &Map, const char *const Label) | |
Epetra_DistObject (const Epetra_DistObject &Source) | |
Epetra_DistObject copy constructor. More... | |
virtual | ~Epetra_DistObject () |
Epetra_DistObject destructor. More... | |
int | Import (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Imports an Epetra_DistObject using the Epetra_Import object. More... | |
int | Import (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Imports an Epetra_DistObject using the Epetra_Export object. More... | |
int | Export (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Exports an Epetra_DistObject using the Epetra_Import object. More... | |
int | Export (const Epetra_SrcDistObject &A, const Epetra_Export &Exporter, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0) |
Exports an Epetra_DistObject using the Epetra_Export object. More... | |
const Epetra_BlockMap & | Map () const |
Returns the address of the Epetra_BlockMap for this multi-vector. More... | |
const Epetra_Comm & | Comm () const |
Returns the address of the Epetra_Comm for this multi-vector. More... | |
bool | DistributedGlobal () const |
Returns true if this multi-vector is distributed global, i.e., not local replicated. More... | |
Public Member Functions inherited from Epetra_Object | |
Epetra_Object (int TracebackModeIn=-1, bool set_label=true) | |
Epetra_Object Constructor. More... | |
Epetra_Object (const char *const Label, int TracebackModeIn=-1) | |
Epetra_Object Constructor. More... | |
Epetra_Object (const Epetra_Object &Object) | |
Epetra_Object Copy Constructor. More... | |
virtual | ~Epetra_Object () |
Epetra_Object Destructor. More... | |
virtual void | SetLabel (const char *const Label) |
Epetra_Object Label definition using char *. More... | |
virtual const char * | Label () const |
Epetra_Object Label access funtion. More... | |
virtual int | ReportError (const std::string Message, int ErrorCode) const |
Error reporting method. More... | |
Public Member Functions inherited from Epetra_SrcDistObject | |
virtual | ~Epetra_SrcDistObject () |
Epetra_SrcDistObject destructor. More... | |
Protected Member Functions | |
int * | All_Indices () const |
int * | IndexOffset () const |
int * | NumIndicesPerRow () const |
int * | NumAllocatedIndicesPerRow () const |
int ** | Indices () const |
int * | Indices (int LocalRow) const |
template<typename int_type > | |
int_type ** | TIndices () const |
template<typename int_type > | |
int_type * | TIndices (int LocalRow) const |
bool | IndicesAreContiguous () const |
bool | StaticProfile () const |
bool | GlobalConstantsComputed () const |
bool | FindGlobalIndexLoc (int LocalRow, int Index, int Start, int &Loc) const |
bool | FindGlobalIndexLoc (int NumIndices, const int *Indices, int Index, int Start, int &Loc) const |
bool | FindGlobalIndexLoc (int LocalRow, long long Index, int Start, int &Loc) const |
bool | FindGlobalIndexLoc (int NumIndices, const long long *Indices, long long Index, int Start, int &Loc) const |
bool | FindMyIndexLoc (int LocalRow, int Index, int Start, int &Loc) const |
bool | FindMyIndexLoc (int NumIndices, const int *Indices, int Index, int Start, int &Loc) const |
int | InsertIndices (int Row, int NumIndices, int *Indices) |
int | InsertIndicesIntoSorted (int Row, int NumIndices, int *Indices) |
int | InsertIndices (int Row, int NumIndices, long long *Indices) |
int | InsertIndicesIntoSorted (int Row, int NumIndices, long long *Indices) |
int | MakeIndicesLocal (const Epetra_BlockMap &DomainMap, const Epetra_BlockMap &RangeMap) |
void | SetIndicesAreLocal (bool Flag) |
void | SetIndicesAreGlobal (bool Flag) |
void | SetSorted (bool Flag) |
int | SortIndices () |
Sort column indices, row-by-row, in ascending order. More... | |
bool | Sorted () const |
If SortIndices() has been called, this query returns true, otherwise it returns false. More... | |
int | RemoveRedundantIndices () |
Removes any redundant column indices in the rows of the graph. More... | |
int | DetermineTriangular () |
bool | NoRedundancies () const |
If RemoveRedundantIndices() has been called, this query returns true, otherwise it returns false. More... | |
Protected Member Functions inherited from Epetra_DistObject | |
virtual int | DoTransfer (const Epetra_SrcDistObject &A, Epetra_CombineMode CombineMode, int NumSameIDs, int NumPermuteIDs, int NumRemoteIDs, int NumExportIDs, int *PermuteToLIDs, int *PermuteFromLIDs, int *RemoteLIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &LenImports, char *&Imports, Epetra_Distributor &Distor, bool DoReverse, const Epetra_OffsetIndex *Indexor) |
Perform actual transfer (redistribution) of data across memory images, using Epetra_Distributor object. More... | |
Protected Member Functions inherited from Epetra_Object | |
std::string | toString (const int &x) const |
std::string | toString (const long long &x) const |
std::string | toString (const double &x) const |
Private Member Functions | |
void | SetGlobalConstantsComputed (bool Flag) |
void | SetIndicesAreContiguous (bool Flag) |
void | SetNoRedundancies (bool Flag) |
void | ComputeIndexState () |
int | MakeColMap (const Epetra_BlockMap &DomainMap, const Epetra_BlockMap &RangeMap) |
int | MakeColMap_int (const Epetra_BlockMap &DomainMap, const Epetra_BlockMap &RangeMap) |
int | MakeColMap_LL (const Epetra_BlockMap &DomainMap, const Epetra_BlockMap &RangeMap) |
int | Allocate (const int *NumIndicesPerRow, int Inc, bool StaticProfile) |
int | ComputeGlobalConstants () |
void | SetFilled (bool Flag) |
bool | Allocated () const |
void | SetAllocated (bool Flag) |
int | CheckSizes (const Epetra_SrcDistObject &A) |
Allows the source and target (this) objects to be compared for compatibility, return nonzero if not. More... | |
int | CopyAndPermute (const Epetra_SrcDistObject &Source, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode=Zero) |
Perform ID copies and permutations that are on processor. More... | |
int | CopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
template<typename int_type > | |
int | CopyAndPermuteRowMatrix (const Epetra_RowMatrix &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
int | CopyAndPermuteCrsGraph (const Epetra_CrsGraph &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
template<typename int_type > | |
int | CopyAndPermuteCrsGraph (const Epetra_CrsGraph &A, int NumSameIDs, int NumPermuteIDs, int *PermuteToLIDs, int *PermuteFromLIDs, const Epetra_OffsetIndex *Indexor, Epetra_CombineMode CombineMode) |
int | PackAndPrepare (const Epetra_SrcDistObject &Source, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor) |
Perform any packing or preparation required for call to DoTransfer(). More... | |
int | PackAndPrepareCrsGraph (const Epetra_CrsGraph &A, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor) |
int | PackAndPrepareRowMatrix (const Epetra_RowMatrix &A, int NumExportIDs, int *ExportLIDs, int &LenExports, char *&Exports, int &SizeOfPacket, int *Sizes, bool &VarSizes, Epetra_Distributor &Distor) |
int | UnpackAndCombine (const Epetra_SrcDistObject &Source, int NumImportIDs, int *ImportLIDs, int LenImports, char *Imports, int &SizeOfPacket, Epetra_Distributor &Distor, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor) |
Perform any unpacking and combining after call to DoTransfer(). More... | |
void | CleanupData () |
template<typename int_type > | |
int | TAllocate (const int *numIndicesPerRow, int Inc, bool staticProfile) |
template<typename int_type > | |
int | InsertGlobalIndices (int_type GlobalRow, int NumIndices, int_type *Indices) |
template<typename int_type > | |
int | TInsertIndices (int Row, int NumIndices, int_type *Indices) |
template<typename int_type > | |
int | TInsertIndicesIntoSorted (int Row, int NumIndices, int_type *Indices) |
template<typename int_type > | |
int | RemoveGlobalIndices (int_type GlobalRow, int NumIndices, int_type *Indices) |
template<typename int_type > | |
int | TRemoveGlobalIndices (long long Row) |
template<typename int_type > | |
bool | FindGlobalIndexLoc (int LocalRow, int_type Index, int Start, int &Loc) const |
template<typename int_type > | |
bool | FindGlobalIndexLoc (int NumIndices, const int_type *Indices, int_type Index, int Start, int &Loc) const |
template<typename int_type > | |
int | ExtractGlobalRowCopy (int_type Row, int LenOfIndices, int &NumIndices, int_type *Indices) const |
template<typename int_type > | |
int | ExtractMyRowCopy (int Row, int LenOfIndices, int &NumIndices, int_type *targIndices) const |
Private Attributes | |
Epetra_CrsGraphData * | CrsGraphData_ |
Friends | |
class | Epetra_CrsMatrix |
class | Epetra_VbrMatrix |
class | Epetra_FECrsGraph |
class | Epetra_FECrsMatrix |
class | Epetra_FEVbrMatrix |
class | Epetra_OffsetIndex |
Insertion/Removal methods | |
int | InsertGlobalIndices (int GlobalRow, int NumIndices, int *Indices) |
Enter a list of elements in a specified global row of the graph. More... | |
int | InsertGlobalIndices (long long GlobalRow, int NumIndices, long long *Indices) |
int | RemoveGlobalIndices (int GlobalRow, int NumIndices, int *Indices) |
Remove a list of elements from a specified global row of the graph. More... | |
int | RemoveGlobalIndices (long long GlobalRow, int NumIndices, long long *Indices) |
int | RemoveGlobalIndices (long long Row) |
Remove all indices from a specified global row of the graph. More... | |
int | InsertMyIndices (int LocalRow, int NumIndices, int *Indices) |
Enter a list of elements in a specified local row of the graph. More... | |
int | RemoveMyIndices (int LocalRow, int NumIndices, int *Indices) |
Remove a list of elements from a specified local row of the graph. More... | |
int | RemoveMyIndices (int Row) |
Remove all indices from a specified local row of the graph. More... | |
Transformation methods | |
int | FillComplete () |
Tranform to local index space. Perform other operations to allow optimal matrix operations. More... | |
int | FillComplete (const Epetra_BlockMap &DomainMap, const Epetra_BlockMap &RangeMap) |
Transform to local index space using specified Domain/Range maps. Perform other operations to allow optimal matrix operations. More... | |
int | OptimizeStorage () |
Make consecutive row index sections contiguous, minimize internal storage used for constructing graph. More... | |
Extraction methods | |
int | ExtractGlobalRowCopy (int GlobalRow, int LenOfIndices, int &NumIndices, int *Indices) const |
Extract a list of elements in a specified global row of the graph. Put into storage allocated by calling routine. More... | |
int | ExtractGlobalRowCopy (long long GlobalRow, int LenOfIndices, int &NumIndices, long long *Indices) const |
int | ExtractMyRowCopy (int LocalRow, int LenOfIndices, int &NumIndices, int *Indices) const |
Extract a list of elements in a specified local row of the graph. Put into storage allocated by calling routine. More... | |
int | ExtractGlobalRowView (int GlobalRow, int &NumIndices, int *&Indices) const |
Get a view of the elements in a specified global row of the graph. More... | |
int | ExtractGlobalRowView (long long GlobalRow, int &NumIndices, long long *&Indices) const |
int | ExtractMyRowView (int LocalRow, int &NumIndices, int *&Indices) const |
Get a view of the elements in a specified local row of the graph. More... | |
Graph Properties Query Methods | |
bool | Filled () const |
If FillComplete() has been called, this query returns true, otherwise it returns false. More... | |
bool | StorageOptimized () const |
If OptimizeStorage() has been called, this query returns true, otherwise it returns false. More... | |
bool | IndicesAreGlobal () const |
If column indices are in global range, this query returns true, otherwise it returns false. More... | |
bool | IndicesAreLocal () const |
If column indices are in local range, this query returns true, otherwise it returns false. More... | |
bool | LowerTriangular () const |
If graph is lower triangular in local index space, this query returns true, otherwise it returns false. More... | |
bool | UpperTriangular () const |
If graph is upper triangular in local index space, this query returns true, otherwise it returns false. More... | |
bool | NoDiagonal () const |
If graph has no diagonal entries in global index space, this query returns true, otherwise it returns false. More... | |
bool | MyGlobalRow (int GID) const |
Returns true of GID is owned by the calling processor, otherwise it returns false. More... | |
bool | MyGlobalRow (long long GID) const |
bool | HaveColMap () const |
Returns true if we have a well-defined ColMap, and returns false otherwise. More... | |
Attribute access functions | |
int | NumMyRows () const |
Returns the number of matrix rows on this processor. More... | |
int | NumGlobalRows () const |
Returns the number of matrix rows in global matrix. More... | |
long long | NumGlobalRows64 () const |
int | NumMyCols () const |
Returns the number of entries in the set of column-indices that appear on this processor. More... | |
int | NumGlobalCols () const |
Returns the number of matrix columns in global matrix. More... | |
long long | NumGlobalCols64 () const |
int | NumGlobalNonzeros () const |
Returns the number of indices in the global graph. More... | |
long long | NumGlobalNonzeros64 () const |
int | NumGlobalDiagonals () const |
Returns the number of diagonal entries in the global graph, based on global row/column index comparisons. More... | |
long long | NumGlobalDiagonals64 () const |
int | NumMyDiagonals () const |
Returns the number of diagonal entries in the local graph, based on global row/column index comparisons. More... | |
int | NumMyBlockRows () const |
Returns the number of block matrix rows on this processor. More... | |
int | NumGlobalBlockRows () const |
Returns the number of Block matrix rows in global matrix. More... | |
long long | NumGlobalBlockRows64 () const |
int | NumMyBlockCols () const |
Returns the number of Block matrix columns on this processor. More... | |
int | NumGlobalBlockCols () const |
Returns the number of Block matrix columns in global matrix. More... | |
long long | NumGlobalBlockCols64 () const |
int | NumMyBlockDiagonals () const |
Returns the number of Block diagonal entries in the local graph, based on global row/column index comparisons. More... | |
int | NumGlobalBlockDiagonals () const |
Returns the number of Block diagonal entries in the global graph, based on global row/column index comparisons. More... | |
long long | NumGlobalBlockDiagonals64 () const |
int | NumGlobalEntries () const |
Returns the number of entries in the global graph. More... | |
long long | NumGlobalEntries64 () const |
int | NumMyEntries () const |
Returns the number of entries on this processor. More... | |
int | MaxRowDim () const |
Returns the max row dimension of block entries on the processor. More... | |
int | GlobalMaxRowDim () const |
Returns the max row dimension of block entries across all processors. More... | |
int | MaxColDim () const |
Returns the max column dimension of block entries on the processor. More... | |
int | GlobalMaxColDim () const |
Returns the max column dimension of block entries across all processors. More... | |
int | NumMyNonzeros () const |
Returns the number of indices in the local graph. More... | |
int | NumGlobalIndices (long long Row) const |
Returns the current number of nonzero entries in specified global row on this processor. More... | |
int | NumAllocatedGlobalIndices (long long Row) const |
Returns the allocated number of nonzero entries in specified global row on this processor. More... | |
int | MaxNumIndices () const |
Returns the maximum number of nonzero entries across all rows on this processor. More... | |
int | GlobalMaxNumIndices () const |
Returns the maximun number of nonzero entries across all rows across all processors. More... | |
int | MaxNumNonzeros () const |
Returns the maximum number of nonzero points across all rows on this processor. More... | |
int | GlobalMaxNumNonzeros () const |
Returns the maximun number of nonzero points across all rows across all processors. More... | |
int | NumMyIndices (int Row) const |
Returns the current number of nonzero entries in specified local row on this processor. More... | |
int | NumAllocatedMyIndices (int Row) const |
Returns the allocated number of nonzero entries in specified local row on this processor. More... | |
int | IndexBase () const |
Returns the index base for row and column indices for this graph. More... | |
long long | IndexBase64 () const |
const Epetra_BlockMap & | RowMap () const |
Returns the RowMap associated with this graph. More... | |
int | ReplaceRowMap (const Epetra_BlockMap &newmap) |
Replaces the current RowMap with the user-specified map object, but only if currentmap->PointSameAs(newmap) is true. More... | |
int | ReplaceColMap (const Epetra_BlockMap &newmap) |
Replaces the current ColMap with the user-specified map object, but only if no entries have been inserted into the graph yet (both IndicesAreLocal() and IndicesAreGlobal() are false) or currentmap->PointSameAs(newmap) is true. More... | |
int | ReplaceDomainMapAndImporter (const Epetra_BlockMap &NewDomainMap, const Epetra_Import *NewImporter) |
Replaces the current DomainMap & Importer with the user-specified map object. More... | |
int | RemoveEmptyProcessesInPlace (const Epetra_BlockMap *NewMap) |
Remove processes owning zero rows from the Maps and their communicator. More... | |
const Epetra_BlockMap & | ColMap () const |
Returns the Column Map associated with this graph. More... | |
const Epetra_BlockMap & | DomainMap () const |
Returns the DomainMap associated with this graph. More... | |
const Epetra_BlockMap & | RangeMap () const |
Returns the RangeMap associated with this graph. More... | |
const Epetra_Import * | Importer () const |
Returns the Importer associated with this graph. More... | |
const Epetra_Export * | Exporter () const |
Returns the Exporter associated with this graph. More... | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this graph. More... | |
Local/Global ID methods | |
int | LRID (int GRID_in) const |
Returns the local row index for given global row index, returns -1 if no local row for this global row. More... | |
int | LRID (long long GRID_in) const |
int | GRID (int LRID_in) const |
Returns the global row index for give local row index, returns IndexBase-1 if we don't have this local row. More... | |
long long | GRID64 (int LRID_in) const |
int | LCID (int GCID_in) const |
Returns the local column index for given global column index, returns -1 if no local column for this global column. More... | |
int | LCID (long long GCID_in) const |
int | GCID (int LCID_in) const |
Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column. More... | |
long long | GCID64 (int LCID_in) const |
bool | MyGRID (int GRID_in) const |
Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
bool | MyGRID (long long GRID_in) const |
bool | MyLRID (int LRID_in) const |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
bool | MyGCID (int GCID_in) const |
Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
bool | MyGCID (long long GCID_in) const |
bool | MyLCID (int LCID_in) const |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false. More... | |
Inlined Operator Methods | |
int * | operator[] (int Loc) |
Inlined bracket operator for fast access to data. (Const and Non-const versions) More... | |
int * | operator[] (int Loc) const |
I/O Methods | |
virtual void | Print (std::ostream &os) const |
Print method. More... | |
void | PrintGraphData (std::ostream &os) const |
void | PrintGraphData (std::ostream &os, int level) const |
Deprecated methods: These methods still work, but will be removed in a future version | |
const Epetra_BlockMap & | ImportMap () const |
Use ColMap() instead. More... | |
int | TransformToLocal () |
Use FillComplete() instead. More... | |
int | TransformToLocal (const Epetra_BlockMap *DomainMap, const Epetra_BlockMap *RangeMap) |
Use FillComplete(const Epetra_BlockMap& DomainMap, const Epetra_BlockMap& RangeMap) instead. More... | |
Expert Users and Developers Only | |
int | ReferenceCount () const |
Returns the reference count of CrsGraphData. More... | |
const Epetra_CrsGraphData * | DataPtr () const |
Returns a pointer to the CrsGraphData instance this CrsGraph uses. More... | |
void | SortGhostsAssociatedWithEachProcessor (bool Flag) |
Forces FillComplete() to locally order ghostnodes associated with each remote processor in ascending order. More... | |
Additional Inherited Members | |
Static Public Member Functions inherited from Epetra_Object | |
static void | SetTracebackMode (int TracebackModeValue) |
Set the value of the Epetra_Object error traceback report mode. More... | |
static int | GetTracebackMode () |
Get the value of the Epetra_Object error report mode. More... | |
static std::ostream & | GetTracebackStream () |
Get the output stream for error reporting. More... | |
Static Public Attributes inherited from Epetra_Object | |
static int | TracebackMode |
Protected Attributes inherited from Epetra_DistObject | |
Epetra_BlockMap | Map_ |
const Epetra_Comm * | Comm_ |
char * | Exports_ |
char * | Imports_ |
int | LenExports_ |
int | LenImports_ |
int * | Sizes_ |
Epetra_CrsGraph: A class for constructing and using sparse compressed row graphs.
Epetra_CrsGraph enables the piecewise construction and use of sparse matrix graphs (the integer structure without values) where entries are intended for row access. Epetra_CrsGraph is an attribute of all Epetra row-based matrix classes, defining their nonzero structure and also holding their Epetra_Map attributes.
Constructing Epetra_CrsGraph objects
Constructing Epetra_CrsGraph objects is a multi-step process. The basic steps are as follows:
Performance Enhancement Issues
The Epetra_CrsGraph class attempts to address four basic types of situations, depending on the user's primary concern:
Notes:
In all but the most advanced uses, users will typically not specify the column map. In other words, graph entries will be submitted using GIDs not LIDs and all entries that are submitted are intended to be inserted into the graph.
If a user is not particularly worried about performance, or really needs the flexibility associated with the first situation, then there is no need to explicitly manage the NumIndicesPerRow values or set StaticProfile to true. In this case, it is best to set NumIndicesPerRow to zero.
Users who are concerned about performance should carefully manage NumIndicesPerRow and set StaticProfile to true. This will give the best performance and use the least amount of memory.
Epetra_Map attributes
Epetra_CrsGraph objects have four Epetra_Map attributes.
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.
Global versus Local indices
After creation and before FillComplete() has been called, the column-indices of the graph are in the global space as received from the user. One of the tasks performed by FillComplete() is to transform the indices to a local index space. The query methods IndicesAreGlobal() and IndicesAreLocal() return true or false depending on whether this transformation has been performed or not.
Note the behavior of several graph methods:
Note that even after a graph is constructed, it is possible to add or remove entries. However, FillComplete must then be called again to restore the graph to a consistent state.
Definition at line 213 of file Epetra_CrsGraph.h.
Epetra_CrsGraph::Epetra_CrsGraph | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | RowMap, | ||
const int * | NumIndicesPerRow, | ||
bool | StaticProfile = false |
||
) |
Epetra_CrsGraph constuctor with variable number of indices per row.
Creates a Epetra_CrsGraph object and allocates storage.
CV | - (In) A Epetra_DataAccess enumerated type set to Copy or View. |
RowMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the rows that this processor will contribute to.In |
NumIndicesPerRow | - (In) An integer array of length NumMyRows such that NumIndicesPerRow[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 66 of file Epetra_CrsGraph.cpp.
Epetra_CrsGraph::Epetra_CrsGraph | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | RowMap, | ||
int | NumIndicesPerRow, | ||
bool | StaticProfile = false |
||
) |
Epetra_CrsGraph constuctor with fixed number of indices per row.
Creates a Epetra_CrsGraph object and allocates storage.
CV | - (In) A Epetra_DataAccess enumerated type set to Copy or View. |
RowMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the rows that this processor will contribute to. |
NumIndicesPerRow | - (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 76 of file Epetra_CrsGraph.cpp.
Epetra_CrsGraph::Epetra_CrsGraph | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | RowMap, | ||
const Epetra_BlockMap & | ColMap, | ||
const int * | NumIndicesPerRow, | ||
bool | StaticProfile = false |
||
) |
Epetra_CrsGraph constuctor with variable number of indices per row.
Creates a Epetra_CrsGraph object and allocates storage.
CV | - (In) A Epetra_DataAccess enumerated type set to Copy or View. |
RowMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the rows that this processor will contribute to. |
ColMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the columns that this processor will contribute to. |
NumIndicesPerRow | - (In) An integer array of length NumMyRows such that NumIndicesPerRow[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 86 of file Epetra_CrsGraph.cpp.
Epetra_CrsGraph::Epetra_CrsGraph | ( | Epetra_DataAccess | CV, |
const Epetra_BlockMap & | RowMap, | ||
const Epetra_BlockMap & | ColMap, | ||
int | NumIndicesPerRow, | ||
bool | StaticProfile = false |
||
) |
Epetra_CrsGraph constuctor with fixed number of indices per row.
Creates a Epetra_CrsGraph object and allocates storage.
CV | - (In) A Epetra_DataAccess enumerated type set to Copy or View. |
RowMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the rows that this processor will contribute to. |
ColMap | - (In) An Epetra_BlockMap (or Epetra_Map or Epetra_LocalMap) listing the columns that this processor will contribute to. |
In | NumIndicesPerRow - 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 100 of file Epetra_CrsGraph.cpp.
Epetra_CrsGraph::Epetra_CrsGraph | ( | const Epetra_CrsGraph & | Graph | ) |
Copy constructor.
This will create a Level 1 deep copy. This Graph will share ownership of the CrsGraphData object with the right hand side Graph.
Definition at line 114 of file Epetra_CrsGraph.cpp.
|
virtual |
Epetra_CrsGraph Destructor.
Definition at line 236 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::InsertGlobalIndices | ( | int | GlobalRow, |
int | NumIndices, | ||
int * | Indices | ||
) |
Enter a list of elements in a specified global row of the graph.
Row | - (In) Global row number of indices. |
NumIndices | - (In) Number of Indices. |
Indices | - (In) Global column indices to insert. |
Definition at line 272 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::InsertGlobalIndices | ( | long long | GlobalRow, |
int | NumIndices, | ||
long long * | Indices | ||
) |
Definition at line 282 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveGlobalIndices | ( | int | GlobalRow, |
int | NumIndices, | ||
int * | Indices | ||
) |
Remove a list of elements from a specified global row of the graph.
Row | - (In) Global row number of indices. |
NumIndices | - (In) Number of Indices. |
Indices | - (In) Global column indices to remove. |
Definition at line 588 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveGlobalIndices | ( | long long | GlobalRow, |
int | NumIndices, | ||
long long * | Indices | ||
) |
Definition at line 598 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveGlobalIndices | ( | long long | Row | ) |
Remove all indices from a specified global row of the graph.
Row | - (In) Global row number of indices. |
Definition at line 701 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::InsertMyIndices | ( | int | LocalRow, |
int | NumIndices, | ||
int * | Indices | ||
) |
Enter a list of elements in a specified local row of the graph.
Row | - (In) Local row number of indices. |
NumIndices | - (In) Number of Indices. |
Indices | - (In) Local column indices to insert. |
Definition at line 291 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveMyIndices | ( | int | LocalRow, |
int | NumIndices, | ||
int * | Indices | ||
) |
Remove a list of elements from a specified local row of the graph.
Row | - (In) Local row number of indices. |
NumIndices | - (In) Number of Indices. |
Indices | - (In) Local column indices to remove. |
Definition at line 607 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveMyIndices | ( | int | Row | ) |
Remove all indices from a specified local row of the graph.
Row | - (In) Local row number of indices. |
Definition at line 721 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::FillComplete | ( | ) |
Tranform to local index space. Perform other operations to allow optimal matrix operations.
This overloading of the FillComplete method assumes that the domain-map and range-map both equal the row-map, and simply calls FillComplete(RowMap(), RowMap()).
Definition at line 963 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::FillComplete | ( | const Epetra_BlockMap & | DomainMap, |
const Epetra_BlockMap & | RangeMap | ||
) |
Transform to local index space using specified Domain/Range maps. Perform other operations to allow optimal matrix operations.
Performs this sequence of operations:
Definition at line 969 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::OptimizeStorage | ( | ) |
Make consecutive row index sections contiguous, minimize internal storage used for constructing graph.
After construction and during initialization (when indices are being added via InsertGlobalIndices() etc.), the column- indices for each row are held in a separate piece of allocated memory. This method moves the column-indices for all rows into one large contiguous array and eliminates internal storage that is not needed after graph 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 indices contiguous and will return a warning code of 1 if the viewed data isn't already contiguous.
Definition at line 1876 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractGlobalRowCopy | ( | int | GlobalRow, |
int | LenOfIndices, | ||
int & | NumIndices, | ||
int * | Indices | ||
) | const |
Extract a list of elements in a specified global row of the graph. Put into storage allocated by calling routine.
Row | - (In) Global row number to get indices. |
LenOfIndices | - (In) Length of Indices array. |
NumIndices | - (Out) Number of Indices. |
Indices | - (Out) Global column indices corresponding to values. |
Definition at line 2049 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractGlobalRowCopy | ( | long long | GlobalRow, |
int | LenOfIndices, | ||
int & | NumIndices, | ||
long long * | Indices | ||
) | const |
Definition at line 2059 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractMyRowCopy | ( | int | LocalRow, |
int | LenOfIndices, | ||
int & | NumIndices, | ||
int * | Indices | ||
) | const |
Extract a list of elements in a specified local row of the graph. Put into storage allocated by calling routine.
Row | - (In) Local row number to get indices. |
LenOfIndices | - (In) Length of Indices array. |
NumIndices | - (Out) Number of Indices. |
Indices | - (Out) Local column indices corresponding to values. |
Definition at line 2091 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractGlobalRowView | ( | int | GlobalRow, |
int & | NumIndices, | ||
int *& | Indices | ||
) | const |
Get a view of the elements in a specified global row of the graph.
This function requires that the graph not be completed (FillComplete() was not called).
Row | - (In) Global row number to get indices. |
NumIndices | - (Out) Number of Indices. |
Indices | - (Out) Global column indices corresponding to values. |
Definition at line 2101 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractGlobalRowView | ( | long long | GlobalRow, |
int & | NumIndices, | ||
long long *& | Indices | ||
) | const |
Definition at line 2123 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ExtractMyRowView | ( | int | LocalRow, |
int & | NumIndices, | ||
int *& | Indices | ||
) | const |
Get a view of the elements in a specified local row of the graph.
This function requires that the graph be completed FillComplete() was called).
Row | - (In) Local row number to get indices. |
NumIndices | - (Out) Number of Indices. |
Indices | - (Out) Column indices corresponding to values. |
Definition at line 2144 of file Epetra_CrsGraph.cpp.
|
inline |
If FillComplete() has been called, this query returns true, otherwise it returns false.
Definition at line 505 of file Epetra_CrsGraph.h.
|
inline |
If OptimizeStorage() has been called, this query returns true, otherwise it returns false.
Definition at line 508 of file Epetra_CrsGraph.h.
|
inline |
If column indices are in global range, this query returns true, otherwise it returns false.
Definition at line 511 of file Epetra_CrsGraph.h.
|
inline |
If column indices are in local range, this query returns true, otherwise it returns false.
Definition at line 514 of file Epetra_CrsGraph.h.
|
inline |
If graph is lower triangular in local index space, this query returns true, otherwise it returns false.
Definition at line 520 of file Epetra_CrsGraph.h.
|
inline |
If graph is upper triangular in local index space, this query returns true, otherwise it returns false.
Definition at line 526 of file Epetra_CrsGraph.h.
|
inline |
If graph has no diagonal entries in global index space, this query returns true, otherwise it returns false.
Definition at line 532 of file Epetra_CrsGraph.h.
|
inline |
Returns true of GID is owned by the calling processor, otherwise it returns false.
Definition at line 536 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 540 of file Epetra_CrsGraph.h.
|
inline |
Returns true if we have a well-defined ColMap, and returns false otherwise.
Definition at line 548 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of matrix rows on this processor.
Definition at line 555 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of matrix rows in global matrix.
Definition at line 559 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 565 of file Epetra_CrsGraph.h.
|
inline |
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.
Definition at line 572 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of matrix columns in global matrix.
Definition at line 579 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 585 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of indices in the global graph.
Note that if the graph's maps are defined such that some nonzeros appear on more than one processor, then those nonzeros will be counted more than once. If the user wishes to assemble a graph from overlapping data, they can use Epetra_FECrsGraph.
Definition at line 595 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 601 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of diagonal entries in the global graph, based on global row/column index comparisons.
Definition at line 608 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 614 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of diagonal entries in the local graph, based on global row/column index comparisons.
Definition at line 620 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of block matrix rows on this processor.
Definition at line 623 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of Block matrix rows in global matrix.
Definition at line 627 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 633 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of Block matrix columns on this processor.
Definition at line 639 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of Block matrix columns in global matrix.
Definition at line 646 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 652 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of Block diagonal entries in the local graph, based on global row/column index comparisons.
Definition at line 658 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of Block diagonal entries in the global graph, based on global row/column index comparisons.
Definition at line 665 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 671 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of entries in the global graph.
Definition at line 678 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 684 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of entries on this processor.
Definition at line 690 of file Epetra_CrsGraph.h.
|
inline |
Returns the max row dimension of block entries on the processor.
Definition at line 695 of file Epetra_CrsGraph.h.
|
inline |
Returns the max row dimension of block entries across all processors.
Definition at line 701 of file Epetra_CrsGraph.h.
|
inline |
Returns the max column dimension of block entries on the processor.
Definition at line 707 of file Epetra_CrsGraph.h.
|
inline |
Returns the max column dimension of block entries across all processors.
Definition at line 713 of file Epetra_CrsGraph.h.
|
inline |
Returns the number of indices in the local graph.
Definition at line 719 of file Epetra_CrsGraph.h.
int Epetra_CrsGraph::NumGlobalIndices | ( | long long | Row | ) | const |
Returns the current number of nonzero entries in specified global row on this processor.
Definition at line 2160 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::NumAllocatedGlobalIndices | ( | long long | Row | ) | const |
Returns the allocated number of nonzero entries in specified global row on this processor.
Definition at line 2173 of file Epetra_CrsGraph.cpp.
|
inline |
Returns the maximum number of nonzero entries across all rows on this processor.
Definition at line 731 of file Epetra_CrsGraph.h.
|
inline |
Returns the maximun number of nonzero entries across all rows across all processors.
Definition at line 737 of file Epetra_CrsGraph.h.
|
inline |
Returns the maximum number of nonzero points across all rows on this processor.
For each entry in the graph, let i = the GRID of the entry and j = the CGID of the entry. Then the entry size is the product of the rowmap elementsize of i and the colmap elementsize of i. Let ki = sum of all entry sizes for the entries in the ith row. For example, if the ith block row had 5 block entries and the element size of each entry was 4-by-4, ki would be 80. Then this function returns the max over all ki for all row on this processor.
Definition at line 749 of file Epetra_CrsGraph.h.
|
inline |
Returns the maximun number of nonzero points across all rows across all processors.
This function returns the max over all processor of MaxNumNonzeros().
Definition at line 756 of file Epetra_CrsGraph.h.
|
inline |
Returns the current number of nonzero entries in specified local row on this processor.
Definition at line 759 of file Epetra_CrsGraph.h.
|
inline |
Returns the allocated number of nonzero entries in specified local row on this processor.
Definition at line 764 of file Epetra_CrsGraph.h.
|
inline |
Returns the index base for row and column indices for this graph.
Index base for this map.
Definition at line 771 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 777 of file Epetra_CrsGraph.h.
|
inline |
Returns the RowMap associated with this graph.
Definition at line 780 of file Epetra_CrsGraph.h.
int Epetra_CrsGraph::ReplaceRowMap | ( | const Epetra_BlockMap & | newmap | ) |
Replaces the current RowMap with the user-specified map object, but only if currentmap->PointSameAs(newmap) is true.
This is a collective function. Returns 0 if map is replaced, -1 if not.
Definition at line 2186 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ReplaceColMap | ( | const Epetra_BlockMap & | newmap | ) |
Replaces the current ColMap with the user-specified map object, but only if no entries have been inserted into the graph yet (both IndicesAreLocal() and IndicesAreGlobal() are false) or currentmap->PointSameAs(newmap) is true.
This is a collective function. Returns 0 if map is replaced, -1 if not.
Definition at line 2199 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::ReplaceDomainMapAndImporter | ( | const Epetra_BlockMap & | NewDomainMap, |
const Epetra_Import * | NewImporter | ||
) |
Replaces the current DomainMap & Importer with the user-specified map object.
Replaces the current DomainMap and Importer with the user-specified map object, but only if the matrix has been FillCompleted, Importer's TargetMap matches the ColMap and Importer's SourceMap matches the DomainMap (assuming the importer isn't null). If an Importer is passed in, Epetra_CrsMatrix will copy it. Returns 0 if map/importer is replaced, -1 if not.
Definition at line 2224 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::RemoveEmptyProcessesInPlace | ( | const Epetra_BlockMap * | NewMap | ) |
Remove processes owning zero rows from the Maps and their communicator.
Remove processes owning zero rows from the Maps and their communicator.
newMap | [in] This must be the result of calling the removeEmptyProcesses() method on the row BlockMap. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes. |
Definition at line 2242 of file Epetra_CrsGraph.cpp.
|
inline |
Returns the Column Map associated with this graph.
Definition at line 830 of file Epetra_CrsGraph.h.
|
inline |
Returns the DomainMap associated with this graph.
Definition at line 836 of file Epetra_CrsGraph.h.
|
inline |
Returns the RangeMap associated with this graph.
Definition at line 842 of file Epetra_CrsGraph.h.
|
inline |
Returns the Importer associated with this graph.
Definition at line 845 of file Epetra_CrsGraph.h.
|
inline |
Returns the Exporter associated with this graph.
Definition at line 848 of file Epetra_CrsGraph.h.
|
inline |
Returns a pointer to the Epetra_Comm communicator associated with this graph.
Definition at line 851 of file Epetra_CrsGraph.h.
|
inline |
Returns the local row index for given global row index, returns -1 if no local row for this global row.
Definition at line 859 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 863 of file Epetra_CrsGraph.h.
|
inline |
Returns the global row index for give local row index, returns IndexBase-1 if we don't have this local row.
Definition at line 874 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 880 of file Epetra_CrsGraph.h.
|
inline |
Returns the local column index for given global column index, returns -1 if no local column for this global column.
Definition at line 887 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 894 of file Epetra_CrsGraph.h.
|
inline |
Returns the global column index for give local column index, returns IndexBase-1 if we don't have this local column.
Definition at line 905 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 911 of file Epetra_CrsGraph.h.
|
inline |
Returns true if the GRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 918 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 922 of file Epetra_CrsGraph.h.
|
inline |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 926 of file Epetra_CrsGraph.h.
|
inline |
Returns true if the GCID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 933 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 937 of file Epetra_CrsGraph.h.
|
inline |
Returns true if the LRID passed in belongs to the calling processor in this map, otherwise returns false.
Definition at line 944 of file Epetra_CrsGraph.h.
|
inline |
Inlined bracket operator for fast access to data. (Const and Non-const versions)
No error checking and dangerous for optimization purposes.
Loc | (In) - Local row. |
Definition at line 957 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 961 of file Epetra_CrsGraph.h.
Epetra_CrsGraph & Epetra_CrsGraph::operator= | ( | const Epetra_CrsGraph & | Source | ) |
Assignment operator.
This will do a Level 1 deep copy. It will share ownership of the CrsGraphData with the right hand side Graph.
Definition at line 3058 of file Epetra_CrsGraph.cpp.
|
virtual |
Print method.
Reimplemented from Epetra_DistObject.
Definition at line 2955 of file Epetra_CrsGraph.cpp.
|
inline |
Definition at line 979 of file Epetra_CrsGraph.h.
|
inline |
Definition at line 980 of file Epetra_CrsGraph.h.
|
inline |
Use ColMap() instead.
Definition at line 987 of file Epetra_CrsGraph.h.
int Epetra_CrsGraph::TransformToLocal | ( | ) |
Use FillComplete() instead.
Definition at line 994 of file Epetra_CrsGraph.cpp.
int Epetra_CrsGraph::TransformToLocal | ( | const Epetra_BlockMap * | DomainMap, |
const Epetra_BlockMap * | RangeMap | ||
) |
Use FillComplete(const Epetra_BlockMap& DomainMap, const Epetra_BlockMap& RangeMap) instead.
Definition at line 999 of file Epetra_CrsGraph.cpp.
|
inline |
Returns the reference count of CrsGraphData.
(Intended for testing purposes.)
Definition at line 1002 of file Epetra_CrsGraph.h.
|
inline |
Returns a pointer to the CrsGraphData instance this CrsGraph uses.
(Intended for developer use only for testing purposes.)
Definition at line 1006 of file Epetra_CrsGraph.h.
|
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.
Definition at line 1018 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1033 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1039 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1042 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1045 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1048 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1051 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1056 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1061 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1067 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1068 of file Epetra_CrsGraph.h.
|
protected |
Definition at line 2897 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 814 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 872 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 826 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 885 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 899 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 932 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 418 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 516 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 429 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 527 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 1767 of file Epetra_CrsGraph.cpp.
|
inlineprotected |
Definition at line 1089 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1090 of file Epetra_CrsGraph.h.
|
inlineprotected |
Definition at line 1091 of file Epetra_CrsGraph.h.
|
protected |
Sort column indices, row-by-row, in ascending order.
Definition at line 1139 of file Epetra_CrsGraph.cpp.
|
inlineprotected |
If SortIndices() has been called, this query returns true, otherwise it returns false.
Definition at line 1101 of file Epetra_CrsGraph.h.
|
protected |
Removes any redundant column indices in the rows of the graph.
Definition at line 1243 of file Epetra_CrsGraph.cpp.
|
protected |
Definition at line 1293 of file Epetra_CrsGraph.cpp.
|
inlineprotected |
If RemoveRedundantIndices() has been called, this query returns true, otherwise it returns false.
Definition at line 1111 of file Epetra_CrsGraph.h.
|
inlineprivate |
Definition at line 1114 of file Epetra_CrsGraph.h.
|
inlineprivate |
Definition at line 1115 of file Epetra_CrsGraph.h.
|
inlineprivate |
Definition at line 1116 of file Epetra_CrsGraph.h.
|
private |
Definition at line 2908 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 1740 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 1348 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 1528 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 202 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 1004 of file Epetra_CrsGraph.cpp.
|
inlineprivate |
Definition at line 1128 of file Epetra_CrsGraph.h.
|
inlineprivate |
Definition at line 1129 of file Epetra_CrsGraph.h.
|
inlineprivate |
Definition at line 1130 of file Epetra_CrsGraph.h.
|
privatevirtual |
Allows the source and target (this) objects to be compared for compatibility, return nonzero if not.
Implements Epetra_DistObject.
Definition at line 2321 of file Epetra_CrsGraph.cpp.
|
privatevirtual |
Perform ID copies and permutations that are on processor.
Implements Epetra_DistObject.
Definition at line 2334 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 2366 of file Epetra_CrsGraph.cpp.
|
private |
|
private |
Definition at line 2466 of file Epetra_CrsGraph.cpp.
|
private |
|
privatevirtual |
Perform any packing or preparation required for call to DoTransfer().
Implements Epetra_DistObject.
Definition at line 2582 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 2652 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 2721 of file Epetra_CrsGraph.cpp.
|
privatevirtual |
Perform any unpacking and combining after call to DoTransfer().
Implements Epetra_DistObject.
Definition at line 2808 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 242 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 123 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 254 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 322 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 441 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 539 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 655 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 760 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 840 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 2018 of file Epetra_CrsGraph.cpp.
|
private |
Definition at line 2070 of file Epetra_CrsGraph.cpp.
|
friend |
Definition at line 1025 of file Epetra_CrsGraph.h.
|
friend |
Definition at line 1026 of file Epetra_CrsGraph.h.
|
friend |
Definition at line 1027 of file Epetra_CrsGraph.h.
|
friend |
Definition at line 1028 of file Epetra_CrsGraph.h.
|
friend |
Definition at line 1029 of file Epetra_CrsGraph.h.
|
friend |
Definition at line 1030 of file Epetra_CrsGraph.h.
|
private |
Definition at line 1215 of file Epetra_CrsGraph.h.