Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
Public Member Functions | Friends | List of all members
Epetra_FECrsGraph Class Reference

#include <Epetra_FECrsGraph.h>

Inheritance diagram for Epetra_FECrsGraph:
Inheritance graph
[legend]
Collaboration diagram for Epetra_FECrsGraph:
Collaboration graph
[legend]

Public Member Functions

 Epetra_FECrsGraph (Epetra_DataAccess CV, const Epetra_BlockMap &RowMap, int *NumIndicesPerRow, bool ignoreNonLocalEntries=false, bool buildNonlocalGraph=false)
 
 Epetra_FECrsGraph (Epetra_DataAccess CV, const Epetra_BlockMap &RowMap, int NumIndicesPerRow, bool ignoreNonLocalEntries=false, bool buildNonlocalGraph=false)
 
 Epetra_FECrsGraph (Epetra_DataAccess CV, const Epetra_BlockMap &RowMap, const Epetra_BlockMap &ColMap, int *NumIndicesPerRow, bool ignoreNonLocalEntries=false, bool buildNonlocalGraph=false)
 
 Epetra_FECrsGraph (Epetra_DataAccess CV, const Epetra_BlockMap &RowMap, const Epetra_BlockMap &ColMap, int NumIndicesPerRow, bool ignoreNonLocalEntries=false, bool buildNonlocalGraph=false)
 
virtual ~Epetra_FECrsGraph ()
 
int InsertGlobalIndices (int numRows, const int *rows, int numCols, const int *cols)
 
int InsertGlobalIndices (int numRows, const long long *rows, int numCols, const long long *cols)
 
int GlobalAssemble (bool callFillComplete=true)
 
int GlobalAssemble (const Epetra_Map &domain_map, const Epetra_Map &range_map, bool callFillComplete=true)
 
bool UseNonlocalGraph () const
 
template<>
std::map< int,
Epetra_CrsGraphData::EntriesInOneRow
< int > > & 
nonlocalRowData ()
 
template<>
std::map< long long,
Epetra_CrsGraphData::EntriesInOneRow
< long long > > & 
nonlocalRowData ()
 
- Public Member Functions inherited from Epetra_CrsGraph
Epetra_CrsGraphoperator= (const Epetra_CrsGraph &Source)
 Assignment operator. More...
 
 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. More...
 
 Epetra_CrsGraph (Epetra_DataAccess CV, const Epetra_BlockMap &RowMap, int NumIndicesPerRow, bool StaticProfile=false)
 Epetra_CrsGraph constuctor with fixed number of indices per row. More...
 
 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. More...
 
 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. More...
 
 Epetra_CrsGraph (const Epetra_CrsGraph &Graph)
 Copy constructor. More...
 
virtual ~Epetra_CrsGraph ()
 Epetra_CrsGraph Destructor.
 
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...
 
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...
 
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...
 
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 column indices are in global range, this query returns true, otherwise it returns false.
 
bool IndicesAreLocal () const
 If column indices are in local range, this query returns true, otherwise it returns false.
 
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.
 
bool MyGlobalRow (long long GID) const
 
bool HaveColMap () const
 Returns true if we have a well-defined ColMap, and returns false otherwise. More...
 
int NumMyRows () const
 Returns the number of matrix rows on this processor.
 
int NumGlobalRows () const
 Returns the number of matrix rows in global matrix.
 
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.
 
int NumGlobalBlockRows () const
 Returns the number of Block matrix rows in global matrix.
 
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.
 
int NumAllocatedGlobalIndices (long long Row) const
 Returns the allocated number of nonzero entries in specified global row on this processor.
 
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.
 
int NumAllocatedMyIndices (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. More...
 
long long IndexBase64 () const
 
const Epetra_BlockMapRowMap () const
 Returns the RowMap associated with this graph.
 
int ReplaceRowMap (const Epetra_BlockMap &newmap)
 
int ReplaceColMap (const Epetra_BlockMap &newmap)
 
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_BlockMapColMap () const
 Returns the Column Map associated with this graph. More...
 
const Epetra_BlockMapDomainMap () const
 Returns the DomainMap associated with this graph. More...
 
const Epetra_BlockMapRangeMap () const
 Returns the RangeMap associated with this graph. More...
 
const Epetra_ImportImporter () const
 Returns the Importer associated with this graph.
 
const Epetra_ExportExporter () const
 Returns the Exporter associated with this graph.
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this graph.
 
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. 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.
 
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. 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...
 
int * operator[] (int Loc)
 Inlined bracket operator for fast access to data. (Const and Non-const versions) More...
 
int * operator[] (int Loc) const
 
virtual void Print (std::ostream &os) const
 Print method.
 
void PrintGraphData (std::ostream &os) const
 
void PrintGraphData (std::ostream &os, int level) const
 
const Epetra_BlockMapImportMap () const
 Use ColMap() instead.
 
int TransformToLocal ()
 Use FillComplete() instead.
 
int TransformToLocal (const Epetra_BlockMap *DomainMap, const Epetra_BlockMap *RangeMap)
 Use FillComplete(const Epetra_BlockMap& DomainMap, const Epetra_BlockMap& RangeMap) instead.
 
int ReferenceCount () const
 Returns the reference count of CrsGraphData. More...
 
const Epetra_CrsGraphDataDataPtr () 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...
 
- 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.
 
virtual ~Epetra_DistObject ()
 Epetra_DistObject destructor.
 
int Import (const Epetra_SrcDistObject &A, const Epetra_Import &Importer, Epetra_CombineMode CombineMode, const Epetra_OffsetIndex *Indexor=0)
 Imports an Epetra_DistObject using the Epetra_Import object. 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_BlockMapMap () const
 Returns the address of the Epetra_BlockMap for this multi-vector.
 
const Epetra_CommComm () const
 Returns the address of the Epetra_Comm for this multi-vector.
 
bool DistributedGlobal () const
 Returns true if this multi-vector is distributed global, i.e., not local replicated.
 
- 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 int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method.
 
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...
 
- Public Member Functions inherited from Epetra_SrcDistObject
virtual ~Epetra_SrcDistObject ()
 Epetra_SrcDistObject destructor.
 

Friends

class Epetra_FECrsMatrix
 

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.
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting.
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 
- Protected Member Functions inherited from Epetra_CrsGraph
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.
 
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.
 
- 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.
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 
- Protected Attributes inherited from Epetra_DistObject
Epetra_BlockMap Map_
 
const Epetra_CommComm_
 
char * Exports_
 
char * Imports_
 
int LenExports_
 
int LenImports_
 
int * Sizes_
 

Detailed Description

Epetra Finite-Element CrsGraph. This class provides the ability to insert indices into a matrix-graph, where the indices represent dense submatrices such as element-stiffnesses that might arise from a finite-element application.

In a parallel setting, indices may be submitted on the local processor for rows that do not reside in the local portion of the row-map. After all indices have been submitted, the GlobalAssemble method gathers all non-local graph rows to the appropriate 'owning' processors (an owning processor is a processor which has the row in its row-map).

Constructor & Destructor Documentation

Epetra_FECrsGraph::Epetra_FECrsGraph ( Epetra_DataAccess  CV,
const Epetra_BlockMap RowMap,
int *  NumIndicesPerRow,
bool  ignoreNonLocalEntries = false,
bool  buildNonlocalGraph = false 
)

Constructor

Epetra_FECrsGraph::Epetra_FECrsGraph ( Epetra_DataAccess  CV,
const Epetra_BlockMap RowMap,
int  NumIndicesPerRow,
bool  ignoreNonLocalEntries = false,
bool  buildNonlocalGraph = false 
)

Constructor

Epetra_FECrsGraph::Epetra_FECrsGraph ( Epetra_DataAccess  CV,
const Epetra_BlockMap RowMap,
const Epetra_BlockMap ColMap,
int *  NumIndicesPerRow,
bool  ignoreNonLocalEntries = false,
bool  buildNonlocalGraph = false 
)

Constructor

Epetra_FECrsGraph::Epetra_FECrsGraph ( Epetra_DataAccess  CV,
const Epetra_BlockMap RowMap,
const Epetra_BlockMap ColMap,
int  NumIndicesPerRow,
bool  ignoreNonLocalEntries = false,
bool  buildNonlocalGraph = false 
)

Constructor

virtual Epetra_FECrsGraph::~Epetra_FECrsGraph ( )
virtual

Constructor Destructor

Member Function Documentation

int Epetra_FECrsGraph::GlobalAssemble ( bool  callFillComplete = true)

Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time. Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method – every processor must enter it before any will complete it.

NOTE***: When GlobalAssemble() calls FillComplete(), it passes the arguments 'DomainMap()' and 'RangeMap()', which are the map attributes held by the base-class CrsMatrix and its graph. If a rectangular matrix is being assembled, the domain-map and range-map must be specified by calling the other overloading of this method. Otherwise, GlobalAssemble() has no way of knowing what these maps should really be.

Parameters
callFillCompleteoption argument, defaults to true. Determines whether GlobalAssemble() internally calls the FillComplete() method on this matrix.
Returns
error-code 0 if successful, non-zero if some error occurs
int Epetra_FECrsGraph::GlobalAssemble ( const Epetra_Map domain_map,
const Epetra_Map range_map,
bool  callFillComplete = true 
)

Gather any overlapping/shared data into the non-overlapping partitioning defined by the Map that was passed to this matrix at construction time. Data imported from other processors is stored on the owning processor with a "sumInto" or accumulate operation. This is a collective method – every processor must enter it before any will complete it.

NOTE***: When GlobalAssemble() (the other overloading of this method) calls FillComplete(), it passes the arguments 'DomainMap()' and 'RangeMap()', which are the map attributes already held by the base-class CrsMatrix and its graph. If a rectangular matrix is being assembled, the domain-map and range-map must be specified. Otherwise, GlobalAssemble() has no way of knowing what these maps should really be.

Parameters
domain_mapuser-supplied domain map for this matrix
range_mapuser-supplied range map for this matrix
callFillCompleteoption argument, defaults to true. Determines whether GlobalAssemble() internally calls the FillComplete() method on this matrix.
Returns
error-code 0 if successful, non-zero if some error occurs
int Epetra_FECrsGraph::InsertGlobalIndices ( int  numRows,
const int *  rows,
int  numCols,
const int *  cols 
)

Insert a rectangular, dense 'submatrix' of entries (matrix nonzero positions) into the graph.

Parameters
numRowsNumber of rows in the submatrix.
rowsList of row-numbers for the submatrix.
numColsNumber of columns in the submatrix.
colsList of column-indices that will be used for each row in the 'rows' list.

The documentation for this class was generated from the following file: