|
Tpetra Matrix/Vector Services
Version of the Day
|
A distributed graph accessed by rows (adjacency lists) and stored sparsely. More...
#include <Tpetra_CrsGraph_decl.hpp>

Public Types | |
| typedef LocalOrdinal | local_ordinal_type |
| This class' first template parameter; the type of local indices. | |
| typedef GlobalOrdinal | global_ordinal_type |
| This class' second template parameter; the type of global indices. | |
| typedef Node | node_type |
| This class' third template parameter; the Kokkos Node type. | |
| typedef Tpetra::Map < LocalOrdinal, GlobalOrdinal, node_type > | map_type |
| The Map specialization used by this class. | |
| typedef Tpetra::Import < LocalOrdinal, GlobalOrdinal, node_type > | import_type |
| The Import specialization used by this class. | |
| typedef Tpetra::Export < LocalOrdinal, GlobalOrdinal, node_type > | export_type |
| The Export specialization used by this class. | |
Typedefs | |
| typedef GlobalOrdinal | packet_type |
| The type of each datum being sent or received in an Import or Export. | |
Public Member Functions | |
Constructor/Destructor Methods | |
| CrsGraph (const Teuchos::RCP< const map_type > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) | |
| Constructor specifying fixed number of entries for each row. | |
| CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) | |
| Constructor specifying (possibly different) number of entries in each row. | |
| CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) | |
| Constructor specifying column Map and fixed number of entries for each row. | |
| CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) | |
| Constructor specifying column Map and number of entries in each row. | |
| CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< LocalOrdinal > &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) | |
| Constructor specifying column Map and arrays containing the graph in sorted local indices. | |
| template<class Node2 > | |
| Teuchos::RCP< CrsGraph < LocalOrdinal, GlobalOrdinal, Node2 > > | clone (const Teuchos::RCP< Node2 > &node2, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) const |
| Create a cloned CrsGraph for a different Node type. | |
| virtual | ~CrsGraph () |
| Destructor. | |
Implementation of Teuchos::ParameterListAcceptor | |
| void | setParameterList (const Teuchos::RCP< Teuchos::ParameterList > ¶ms) |
| Set the given list of parameters (must be nonnull). | |
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () const |
| Default parameter list suitable for validation. | |
Insertion/Removal Methods | |
| void | insertGlobalIndices (GlobalOrdinal globalRow, const Teuchos::ArrayView< const GlobalOrdinal > &indices) |
| Insert global indices into the graph. | |
| void | insertLocalIndices (const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &indices) |
| Insert local indices into the graph. | |
| void | removeLocalIndices (LocalOrdinal localRow) |
| Remove all graph indices from the specified local row. | |
Transformational Methods | |
| void | globalAssemble () |
| Communicate non-local contributions to other processes. | |
| void | resumeFill (const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) |
| void | fillComplete (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) |
| Signal that data entry is complete, specifying domain and range maps. | |
| void | fillComplete (const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null) |
| Signal that data entry is complete. | |
| void | expertStaticFillComplete (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer=Teuchos::null, const Teuchos::RCP< const export_type > &exporter=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null) |
| Perform a fillComplete on a graph that already has data, via setAllIndices(). | |
Methods implementing RowGraph. | |
| Teuchos::RCP< const Comm< int > > | getComm () const |
| Returns the communicator. | |
| Teuchos::RCP< Node > | getNode () const |
| Returns the underlying node. | |
| Teuchos::RCP< const map_type > | getRowMap () const |
| Returns the Map that describes the row distribution in this graph. | |
| Teuchos::RCP< const map_type > | getColMap () const |
| Returns the Map that describes the column distribution in this graph. | |
| Teuchos::RCP< const map_type > | getDomainMap () const |
| Returns the Map associated with the domain of this graph. | |
| Teuchos::RCP< const map_type > | getRangeMap () const |
| Returns the Map associated with the domain of this graph. | |
| Teuchos::RCP< const import_type > | getImporter () const |
| Returns the importer associated with this graph. | |
| Teuchos::RCP< const export_type > | getExporter () const |
| Returns the exporter associated with this graph. | |
| global_size_t | getGlobalNumRows () const |
| Returns the number of global rows in the graph. | |
| global_size_t | getGlobalNumCols () const |
| Returns the number of global columns in the graph. | |
| size_t | getNodeNumRows () const |
| Returns the number of graph rows owned on the calling node. | |
| size_t | getNodeNumCols () const |
| Returns the number of columns connected to the locally owned rows of this graph. | |
| GlobalOrdinal | getIndexBase () const |
| Returns the index base for global indices for this graph. | |
| global_size_t | getGlobalNumEntries () const |
| Returns the global number of entries in the graph. | |
| size_t | getNodeNumEntries () const |
| Returns the local number of entries in the graph. | |
| size_t | getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const |
| Returns the current number of entries on this node in the specified global row. | |
| size_t | getNumEntriesInLocalRow (LocalOrdinal localRow) const |
| Returns the current number of entries on this node in the specified local row. | |
| size_t | getNodeAllocationSize () const |
| Returns the total number of indices allocated for the graph, across all rows on this node. | |
| size_t | getNumAllocatedEntriesInGlobalRow (GlobalOrdinal globalRow) const |
| Returns the current number of allocated entries for this node in the specified global row . | |
| size_t | getNumAllocatedEntriesInLocalRow (LocalOrdinal localRow) const |
| Returns the current number of allocated entries on this node in the specified local row. | |
| global_size_t | getGlobalNumDiags () const |
| Returns the number of global diagonal entries, based on global row/column index comparisons. | |
| size_t | getNodeNumDiags () const |
| Returns the number of local diagonal entries, based on global row/column index comparisons. | |
| size_t | getGlobalMaxNumRowEntries () const |
| Maximum number of entries in all rows over all processes. | |
| size_t | getNodeMaxNumRowEntries () const |
| Maximum number of entries in all rows owned by the calling process. | |
| bool | hasColMap () const |
| Whether the graph has a column Map. | |
| bool | isLowerTriangular () const |
| Whether the graph is locally lower triangular. | |
| bool | isUpperTriangular () const |
| Whether the graph is locally upper triangular. | |
| bool | isLocallyIndexed () const |
| Whether column indices are stored using local indices on the calling process. | |
| bool | isGloballyIndexed () const |
| Whether column indices are stored using global indices on the calling process. | |
| bool | isFillComplete () const |
| Whether fillComplete() has been called and the graph is in compute mode. | |
| bool | isFillActive () const |
| Whether resumeFill() has been called and the graph is in edit mode. | |
| bool | isSorted () const |
| Whether graph indices in all rows are known to be sorted. | |
| bool | isStorageOptimized () const |
Returns true if storage has been optimized. | |
| ProfileType | getProfileType () const |
Returns true if the graph was allocated with static data structures. | |
| void | getGlobalRowCopy (GlobalOrdinal GlobalRow, const Teuchos::ArrayView< GlobalOrdinal > &Indices, size_t &NumIndices) const |
| Get a copy of the column indices (as global indices) in the given row. | |
| void | getLocalRowCopy (LocalOrdinal LocalRow, const Teuchos::ArrayView< LocalOrdinal > &indices, size_t &NumIndices) const |
| Get a copy of the column indices (as local indices) in the given row. | |
| void | getGlobalRowView (GlobalOrdinal GlobalRow, Teuchos::ArrayView< const GlobalOrdinal > &Indices) const |
| Return a const, nonpersisting view of global indices in the given row. | |
| void | getLocalRowView (LocalOrdinal LocalRow, Teuchos::ArrayView< const LocalOrdinal > &indices) const |
| Return a const, nonpersisting view of local indices in the given row. | |
Overridden from Teuchos::Describable | |
| std::string | description () const |
| Return a simple one-line description of this object. | |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const |
| Print the object with some verbosity level to an FancyOStream object. | |
Implementation of DistObject | |
| virtual bool | checkSizes (const SrcDistObject &source) |
| Compare the source and target (this) objects for compatibility. | |
| virtual void | copyAndPermute (const SrcDistObject &source, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs) |
| Perform copies and permutations that are local to this process. | |
| virtual void | packAndPrepare (const SrcDistObject &source, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< GlobalOrdinal > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) |
| Perform any packing or preparation required for communication. | |
| virtual void | pack (const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Teuchos::Array< GlobalOrdinal > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, Distributor &distor) const |
| Pack this object's data for Import or Export. | |
| virtual void | unpackAndCombine (const Teuchos::ArrayView< const LocalOrdinal > &importLIDs, const Teuchos::ArrayView< const GlobalOrdinal > &imports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t constantNumPackets, Distributor &distor, CombineMode CM) |
| Perform any unpacking and combining after communication. | |
Advanced methods, at increased risk of deprecation. | |
| void | getNumEntriesPerLocalRowUpperBound (Teuchos::ArrayRCP< const size_t > &boundPerLocalRow, size_t &boundForAllLocalRows, bool &boundSameForAllLocalRows) const |
| Get an upper bound on the number of entries that can be stored in each row. | |
| void | setAllIndices (const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< LocalOrdinal > &columnIndices) |
| Set the graph's data directly, using 1-D storage. | |
| Teuchos::ArrayRCP< const size_t > | getNodeRowPtrs () const |
| Get an ArrayRCP of the row-offsets. | |
| Teuchos::ArrayRCP< const LocalOrdinal > | getNodePackedIndices () const |
| Get an ArrayRCP of the packed column-indices. | |
| void | replaceColMap (const Teuchos::RCP< const map_type > &newColMap) |
| Replace the graph's current column Map with the given Map. | |
| void | reindexColumns (const Teuchos::RCP< const map_type > &newColMap, const Teuchos::RCP< const import_type > &newImport=Teuchos::null, const bool sortIndicesInEachRow=true) |
| Reindex the column indices in place, and replace the column Map. Optionally, replace the Import object as well. | |
| void | replaceDomainMapAndImporter (const Teuchos::RCP< const map_type > &newDomainMap, const Teuchos::RCP< const import_type > &newImporter) |
| Replace the current domain Map and Import with the given parameters. | |
| virtual void | removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) |
| Remove processes owning zero rows from the Maps and their communicator. | |
Public methods for redistributing data | |
| void | doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM) |
| Import data into this object using an Import object ("forward mode"). | |
| void | doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM) |
| Import data into this object using an Export object ("reverse mode"). | |
| void | doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM) |
| Export data into this object using an Export object ("forward mode"). | |
| void | doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM) |
| Export data into this object using an Import object ("reverse mode"). | |
Attribute accessor methods | |
| bool | isDistributed () const |
| Whether this is a globally distributed object. | |
| virtual Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > | getMap () const |
| The Map describing the parallel distribution of this object. | |
I/O methods | |
| void | print (std::ostream &os) const |
| Print this object to the given output stream. | |
Protected Member Functions | |
| void | setDomainRangeMaps (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap) |
| void | checkInternalState () const |
| Check whether the graph's state is valid. | |
| bool | hasRowInfo () const |
| Whether it is valid to call getRowInfo(). | |
| virtual size_t | constantNumberOfPackets () const |
| Whether the implementation's instance promises always to have a constant number of packets per LID, and if so, how many packets per LID there are. | |
| virtual void | doTransfer (const SrcDistObject &src, CombineMode CM, size_t numSameIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteToLIDs, const Teuchos::ArrayView< const LocalOrdinal > &permuteFromLIDs, const Teuchos::ArrayView< const LocalOrdinal > &remoteLIDs, const Teuchos::ArrayView< const LocalOrdinal > &exportLIDs, Distributor &distor, ReverseOption revOp) |
| Redistribute data across memory images. | |
| virtual void | createViews () const |
| Hook for creating a const view. | |
| virtual void | createViewsNonConst (KokkosClassic::ReadWriteOption rwo) |
| Hook for creating a nonconst view. | |
| virtual void | releaseViews () const |
| Hook for releasing views. | |
Methods governing changes between global and local indices | |
| void | makeColMap () |
| Make the graph's column Map, if it does not already have one. | |
| void | makeIndicesLocal () |
| void | makeImportExport () |
Methods for inserting indices or transforming values | |
| template<ELocalGlobal lg> | |
| size_t | filterIndices (const SLocalGlobalNCViews &inds) const |
| template<class T > | |
| size_t | filterGlobalIndicesAndValues (const Teuchos::ArrayView< GlobalOrdinal > &ginds, const Teuchos::ArrayView< T > &vals) const |
| template<class T > | |
| size_t | filterLocalIndicesAndValues (const Teuchos::ArrayView< LocalOrdinal > &linds, const Teuchos::ArrayView< T > &vals) const |
| size_t | insertIndices (const RowInfo &rowInfo, const SLocalGlobalViews &newInds, const ELocalGlobal lg, const ELocalGlobal I) |
| Insert indices into the given row. | |
| template<class Scalar > | |
| void | insertIndicesAndValues (const RowInfo &rowInfo, const SLocalGlobalViews &newInds, const Teuchos::ArrayView< Scalar > &oldRowVals, const Teuchos::ArrayView< const Scalar > &newRowVals, const ELocalGlobal lg, const ELocalGlobal I) |
| Insert indices and their values into the given row. | |
| void | insertGlobalIndicesImpl (const LocalOrdinal myRow, const Teuchos::ArrayView< const GlobalOrdinal > &indices) |
| void | insertLocalIndicesImpl (const LocalOrdinal myRow, const Teuchos::ArrayView< const LocalOrdinal > &indices) |
| void | insertLocalIndicesFiltered (const LocalOrdinal localRow, const Teuchos::ArrayView< const LocalOrdinal > &indices) |
| Like insertLocalIndices(), but with column Map filtering. | |
| void | insertGlobalIndicesFiltered (const GlobalOrdinal localRow, const Teuchos::ArrayView< const GlobalOrdinal > &indices) |
| Like insertGlobalIndices(), but with column Map filtering. | |
| template<class Scalar , class BinaryFunction > | |
| LocalOrdinal | transformLocalValues (RowInfo rowInfo, const Teuchos::ArrayView< Scalar > &rowVals, const Teuchos::ArrayView< const LocalOrdinal > &inds, const Teuchos::ArrayView< const Scalar > &newVals, BinaryFunction f) const |
| Transform the given values using local indices. | |
| template<class Scalar , class BinaryFunction > | |
| LocalOrdinal | transformGlobalValues (RowInfo rowInfo, const Teuchos::ArrayView< Scalar > &rowVals, const Teuchos::ArrayView< const GlobalOrdinal > &inds, const Teuchos::ArrayView< const Scalar > &newVals, BinaryFunction f) const |
| Transform the given values using global indices. | |
Methods for sorting and merging column indices. | |
| bool | isMerged () const |
| Whether duplicate column indices in each row have been merged. | |
| void | setLocallyModified () |
| Report that we made a local modification to its structure. | |
| void | sortAllIndices () |
| Sort the column indices in all the rows. | |
| void | sortRowIndices (RowInfo rowinfo) |
| Sort the column indices in the given row. | |
| template<class Scalar > | |
| void | sortRowIndicesAndValues (const RowInfo rowinfo, const Teuchos::ArrayView< Scalar > &values) |
| Sort the column indices and their values in the given row. | |
| void | mergeAllIndices () |
| Merge duplicate row indices in all of the rows. | |
| void | mergeRowIndices (RowInfo rowinfo) |
| Merge duplicate row indices in the given row. | |
| template<class Scalar > | |
| void | mergeRowIndicesAndValues (RowInfo rowinfo, const Teuchos::ArrayView< Scalar > &rowValues) |
| Merge duplicate row indices in the given row, along with their corresponding values. | |
Methods to access the graph's data | |
| RowInfo | getRowInfo (size_t myRow) const |
| Return indexing information for the given (local) row index. | |
| Teuchos::ArrayView< const LocalOrdinal > | getLocalView (RowInfo rowinfo) const |
| Get a const view of the local column indices in the given row. | |
| Teuchos::ArrayView< LocalOrdinal > | getLocalViewNonConst (RowInfo rowinfo) |
| Get a nonconst view of the local column indices in the given row. | |
| Teuchos::ArrayView< const GlobalOrdinal > | getGlobalView (RowInfo rowinfo) const |
| Get a const view of the global column indices in the given row. | |
| Teuchos::ArrayView< GlobalOrdinal > | getGlobalViewNonConst (RowInfo rowinfo) |
| Get a nonconst view of the global column indices in the given row. | |
| size_t | findLocalIndex (RowInfo rowinfo, LocalOrdinal ind, size_t hint=0) const |
| Find the column offset corresponding to the given (local) column index. | |
| size_t | findLocalIndex (RowInfo rowinfo, LocalOrdinal ind, Teuchos::ArrayView< const LocalOrdinal > colInds, size_t hint=0) const |
| size_t | findGlobalIndex (RowInfo rowinfo, GlobalOrdinal ind, size_t hint=0) const |
| Find the column offset corresponding to the given (global) column index. | |
Methods to set or access the local Kokkos graph. | |
| void | fillLocalGraph (const Teuchos::RCP< Teuchos::ParameterList > ¶ms) |
|
const Teuchos::RCP< const local_graph_type > | getLocalGraph () const |
|
const Teuchos::RCP < local_graph_type > | getLocalGraphNonConst () |
Protected Attributes | |
| Teuchos::RCP< const map_type > | rowMap_ |
| The Map describing the distribution of rows of the graph. | |
| Teuchos::RCP< const map_type > | colMap_ |
| The Map describing the distribution of columns of the graph. | |
| Teuchos::RCP< const map_type > | rangeMap_ |
| The Map describing the range of the (matrix corresponding to the) graph. | |
| Teuchos::RCP< const map_type > | domainMap_ |
| The Map describing the domain of the (matrix corresponding to the) graph. | |
| Teuchos::RCP< const import_type > | importer_ |
| The Import from the domain Map to the column Map. | |
| Teuchos::RCP< const export_type > | exporter_ |
| The Export from the row Map to the range Map. | |
| ProfileType | pftype_ |
| Whether the graph was allocated with static or dynamic profile. | |
| Teuchos::ArrayRCP< const size_t > | numAllocPerRow_ |
| The maximum number of entries to allow in each locally owned row, per row. | |
| size_t | numAllocForAllRows_ |
| The maximum number of entries to allow in each locally owned row. | |
| Teuchos::ArrayRCP< LocalOrdinal > | lclInds1D_ |
| lclInds1D_ are the indices for all rows | |
| Teuchos::ArrayRCP< GlobalOrdinal > | gblInds1D_ |
| gblInds1D_ are the indices for all rows | |
| Teuchos::ArrayRCP < Teuchos::Array< LocalOrdinal > > | lclInds2D_ |
lclInds2D_[r] are the indices for row r. | |
| Teuchos::ArrayRCP < Teuchos::Array < GlobalOrdinal > > | gblInds2D_ |
gblInds2D_[r] are the indices for row r. | |
| Teuchos::ArrayRCP< size_t > | numRowEntries_ |
| The number of local entries in each locally owned row. | |
| bool | lowerTriangular_ |
| Whether the graph is locally lower triangular. | |
| bool | upperTriangular_ |
| Whether the graph is locally upper triangular. | |
| bool | indicesAreSorted_ |
| Whether the graph's indices are sorted in each row, on this process. | |
| bool | noRedundancies_ |
| Whether the graph's indices are non-redundant (merged) in each row, on this process. | |
| bool | haveLocalConstants_ |
| Whether this process has computed local constants. | |
| bool | haveGlobalConstants_ |
| Whether all processes have computed global constants. | |
| std::map< GlobalOrdinal, std::vector< GlobalOrdinal > > | nonlocals_ |
| Nonlocal data given to insertGlobalValues() or sumIntoGlobalValues(). | |
| bool | sortGhostsAssociatedWithEachProcessor_ |
| Whether to require makeColMap() (and therefore fillComplete()) to order column Map GIDs associated with each remote process in ascending order. | |
| Teuchos::RCP< const Map < LocalOrdinal, GlobalOrdinal, Node > > | map_ |
| The Map over which this object is distributed. | |
| Teuchos::Array< GlobalOrdinal > | imports_ |
| Buffer into which packed data are imported (received from other processes). | |
| Teuchos::Array< size_t > | numImportPacketsPerLID_ |
| Number of packets to receive for each receive operation. | |
| Teuchos::Array< GlobalOrdinal > | exports_ |
| Buffer from which packed data are exported (sent to other processes). | |
| Teuchos::Array< size_t > | numExportPacketsPerLID_ |
| Number of packets to send for each send operation. | |
Related Functions | |
(Note that these are not member functions.) | |
| template<class LocalOrdinal , class GlobalOrdinal , class Node > | |
| Teuchos::RCP< CrsGraph < LocalOrdinal, GlobalOrdinal, Node > > | createCrsGraph (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null) |
| Non-member function to create an empty CrsGraph given a row map and a non-zero profile. | |
Methods for use only by experts | |
| void | removeEmptyProcessesInPlace (Teuchos::RCP< Tpetra::DistObject< PT, LO, GO, NT > > &input, const Teuchos::RCP< const Map< LO, GO, NT > > &newMap) |
| void | removeEmptyProcessesInPlace (Teuchos::RCP< Tpetra::DistObject< PT, LO, GO, NT > > &input) |
| virtual void | removeEmptyProcessesInPlace (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap) |
| Remove processes which contain no elements in this object's Map. | |
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
| LocalOrdinal | The type of local indices. See the documentation of Map for requirements. |
| GlobalOrdinal | The type of global indices. See the documentation of Map for requirements. |
| Node | The Kokkos Node type. See the documentation of Map for requirements. |
This class implements a distributed-memory parallel sparse graph. It provides access by rows to the elements of the graph, as if the local data were stored in compressed sparse row format (adjacency lists, in graph terms). (Implementations are not required to store the data in this way internally.) This class has an interface like that of Epetra_CrsGraph, but also allows insertion of data into nonowned rows, much like Epetra_FECrsGraph.
Before reading the rest of this documentation, it helps to know something about the Teuchos memory management classes, in particular Teuchos::RCP, Teuchos::ArrayRCP, and Teuchos::ArrayView. You should also know a little bit about MPI (the Message Passing Interface for distributed-memory programming). You won't have to use MPI directly to use CrsGraph, but it helps to be familiar with the general idea of distributed storage of data over a communicator. Finally, you should read the documentation of Map.
Graph entries can be added using either local or global coordinates for the indices. The accessors isGloballyIndexed() and isLocallyIndexed() indicate whether the indices are currently stored as global or local indices. Many of the class methods are divided into global and local versions, which differ only in whether they accept/return indices in the global or local coordinate space. Some of these methods may only be used if the graph coordinates are in the appropriate coordinates. For example, getGlobalRowView() returns a View to the indices in global coordinates; if the indices are not in global coordinates, then no such View can be created.
The global/local distinction does distinguish between operation on the global/local graph. Almost all methods operate on the local graph, i.e., the rows of the graph associated with the local node, per the distribution specified by the row map. Access to non-local rows requires performing an explicit communication via the import/export capabilities of the CrsGraph object; see DistObject. However, the method insertGlobalIndices() is an exception to this rule, as non-local rows are allowed to be added via the local graph. These rows are stored in the local graph and communicated to the appropriate node on the next call to globalAssemble() or fillComplete() (the latter calls the former).
Definition at line 173 of file Tpetra_CrsGraph_decl.hpp.
| typedef LocalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type |
This class' first template parameter; the type of local indices.
Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 189 of file Tpetra_CrsGraph_decl.hpp.
| typedef GlobalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type |
This class' second template parameter; the type of global indices.
Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 191 of file Tpetra_CrsGraph_decl.hpp.
| typedef Node Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::node_type |
This class' third template parameter; the Kokkos Node type.
Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 193 of file Tpetra_CrsGraph_decl.hpp.
| typedef Tpetra::Map<LocalOrdinal, GlobalOrdinal, node_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type |
The Map specialization used by this class.
Definition at line 196 of file Tpetra_CrsGraph_decl.hpp.
| typedef Tpetra::Import<LocalOrdinal, GlobalOrdinal, node_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::import_type |
The Import specialization used by this class.
Definition at line 198 of file Tpetra_CrsGraph_decl.hpp.
| typedef Tpetra::Export<LocalOrdinal, GlobalOrdinal, node_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::export_type |
The Export specialization used by this class.
Definition at line 200 of file Tpetra_CrsGraph_decl.hpp.
typedef GlobalOrdinal Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::packet_type [inherited] |
The type of each datum being sent or received in an Import or Export.
Note that this type does not always correspond to the Scalar template parameter of subclasses.
Definition at line 183 of file Tpetra_DistObject_decl.hpp.
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph | ( | const Teuchos::RCP< const map_type > & | rowMap, |
| size_t | maxNumEntriesPerRow, | ||
| ProfileType | pftype = DynamicProfile, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Constructor specifying fixed number of entries for each row.
| rowMap | [in] Distribution of rows of the graph. |
| maxNumEntriesPerRow | [in] Maximum number of graph entries per row. If pftype==DynamicProfile, this is only a hint, and you can set this to zero without affecting correctness. If pftype==StaticProfile, this sets the amount of storage allocated, and you cannot exceed this number of entries in any row. |
| pftype | [in] Whether to allocate storage dynamically (DynamicProfile) or statically (StaticProfile). |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. |
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph | ( | const Teuchos::RCP< const map_type > & | rowMap, |
| const Teuchos::ArrayRCP< const size_t > & | NumEntriesPerRowToAlloc, | ||
| ProfileType | pftype = DynamicProfile, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Constructor specifying (possibly different) number of entries in each row.
| rowMap | [in] Distribution of rows of the graph. |
| NumEntriesPerRowToAlloc | [in] Maximum number of graph entries to allocate for each row. If pftype==DynamicProfile, this is only a hint. If pftype==StaticProfile, this sets the amount of storage allocated, and you cannot exceed the allocated number of entries for any row. |
| pftype | [in] Whether to allocate storage dynamically (DynamicProfile) or statically (StaticProfile). |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. |
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph | ( | const Teuchos::RCP< const map_type > & | rowMap, |
| const Teuchos::RCP< const map_type > & | colMap, | ||
| size_t | maxNumEntriesPerRow, | ||
| ProfileType | pftype = DynamicProfile, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Constructor specifying column Map and fixed number of entries for each row.
| rowMap | [in] Distribution of rows of the graph. |
| colMap | [in] Distribution of columns of the graph. |
| maxNumEntriesPerRow | [in] Maximum number of graph entries per row. If pftype==DynamicProfile, this is only a hint, and you can set this to zero without affecting correctness. If pftype==StaticProfile, this sets the amount of storage allocated, and you cannot exceed this number of entries in any row. |
| pftype | [in] Whether to allocate storage dynamically (DynamicProfile) or statically (StaticProfile). |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. |
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph | ( | const Teuchos::RCP< const map_type > & | rowMap, |
| const Teuchos::RCP< const map_type > & | colMap, | ||
| const Teuchos::ArrayRCP< const size_t > & | NumEntriesPerRowToAlloc, | ||
| ProfileType | pftype = DynamicProfile, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Constructor specifying column Map and number of entries in each row.
| rowMap | [in] Distribution of rows of the graph. |
| colMap | [in] Distribution of columns of the graph. |
| NumEntriesPerRowToAlloc | [in] Maximum number of graph entries to allocate for each row. If pftype==DynamicProfile, this is only a hint. If pftype==StaticProfile, this sets the amount of storage allocated, and you cannot exceed the allocated number of entries for any row. |
| pftype | [in] Whether to allocate storage dynamically (DynamicProfile) or statically (StaticProfile). |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. |
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph | ( | const Teuchos::RCP< const map_type > & | rowMap, |
| const Teuchos::RCP< const map_type > & | colMap, | ||
| const Teuchos::ArrayRCP< size_t > & | rowPointers, | ||
| const Teuchos::ArrayRCP< LocalOrdinal > & | columnIndices, | ||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Constructor specifying column Map and arrays containing the graph in sorted local indices.
| rowMap | [in] Distribution of rows of the graph. |
| colMap | [in] Distribution of columns of the graph. |
| rowPointers | [in] The beginning of each row in the graph, as in a CSR "rowptr" array. The length of this vector should be equal to the number of rows in the graph, plus one. This last entry should store the nunber of nonzeros in the graph. |
| columnIndices | [in] The local indices of the columns, as in a CSR "colind" array. The length of this vector should be equal to the number of unknowns in the graph. |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. |
| Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::~CrsGraph | ( | ) | [virtual] |
Destructor.
Definition at line 251 of file Tpetra_CrsGraph_def.hpp.
| Teuchos::RCP<CrsGraph<LocalOrdinal, GlobalOrdinal, Node2> > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::clone | ( | const Teuchos::RCP< Node2 > & | node2, |
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) | const [inline] |
Create a cloned CrsGraph for a different Node type.
This method creates a new CrsGraph on a specified Kokkos Node type, with all of the entries of this CrsGraph object.
| node2 | [in] Kokkos Node instance for constructing the clone CrsGraph and its constituent objects. |
| params | [in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values. See the list below for valid options. |
Parameters accepted by this method:
true, creates the clone with a static allocation profile. If false, a dynamic allocation profile is used.true, fills clone using this graph's column map and local indices (requires that this graph have a column map.) If false, fills clone using global indices and does not provide a column map. By default, will use local indices only if this graph is using local indices.true, calls fillComplete() on the cloned CrsGraph object, with parameters from params sublist "CrsGraph". The domain map and range maps passed to fillComplete() are those of the map being cloned, if they exist. Otherwise, the row map is used. Definition at line 354 of file Tpetra_CrsGraph_decl.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setParameterList | ( | const Teuchos::RCP< Teuchos::ParameterList > & | params | ) |
Set the given list of parameters (must be nonnull).
Definition at line 289 of file Tpetra_CrsGraph_def.hpp.
| RCP< const ParameterList > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getValidParameters | ( | ) | const |
Default parameter list suitable for validation.
Definition at line 257 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndices | ( | GlobalOrdinal | globalRow, |
| const Teuchos::ArrayView< const GlobalOrdinal > & | indices | ||
| ) |
Insert global indices into the graph.
globalRow is a valid index in the row Map. It need not be owned by the calling process. isLocallyIndexed() == false isStorageOptimized() == falseindicesAreAllocated() == true isGloballyIndexed() == trueIf globalRow does not belong to the graph on this process, then it will be communicated to the appropriate process when globalAssemble() is called. (That method will be called automatically during the next call to fillComplete().) Otherwise, the entries will be inserted into the part of the graph owned by the calling process.
If the graph row already contains entries at the indices corresponding to values in indices, then the redundant indices will be eliminated. This may happen either at insertion or during the next call to fillComplete().
Definition at line 2027 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertLocalIndices | ( | const LocalOrdinal | localRow, |
| const Teuchos::ArrayView< const LocalOrdinal > & | indices | ||
| ) |
Insert local indices into the graph.
localRow is a local row belonging to the graph on this process. isGloballyIndexed() == false isStorageOptimized() == false hasColMap() == trueindicesAreAllocated() == true isLocallyIndexed() == trueindices, then the redundant indices will be eliminated; this may happen at insertion or during the next call to fillComplete(). Definition at line 1898 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::removeLocalIndices | ( | LocalOrdinal | localRow | ) |
Remove all graph indices from the specified local row.
localRow is a local row of this graph. isGloballyIndexed() == false isStorageOptimized() == falsegetNumEntriesInLocalRow(localRow) == 0 indicesAreAllocated() == true isLocallyIndexed() == true Definition at line 2170 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::globalAssemble | ( | ) |
Communicate non-local contributions to other processes.
Each of the methods in this group is a global collective. It is necessary to call these mehtods on all nodes participating in the communicator associated with this graph. This method is called automatically by fillComplete(). Most users do not need to call this themselves, though we do permit this.
Definition at line 2311 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::resumeFill | ( | const Teuchos::RCP< Teuchos::ParameterList > & | params = null | ) |
Resume fill operations. After calling fillComplete(), resumeFill() must be called before initiating any changes to the graph.
resumeFill() may be called repeatedly.
isFillActive() == true isFillComplete() == false Definition at line 2597 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::fillComplete | ( | const Teuchos::RCP< const map_type > & | domainMap, |
| const Teuchos::RCP< const map_type > & | rangeMap, | ||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = null |
||
| ) |
Signal that data entry is complete, specifying domain and range maps.
Off-process indices are distributed (via globalAssemble()), indices are sorted, redundant indices are eliminated, and global indices are transformed to local indices.
isFillActive() == true isFillComplete()() == false isFillActive() == false isFillComplete() == trueParameters:
"Optimize Storage" (bool): Default is false. If true, then isStorageOptimized() returns true after fillComplete finishes. See isStorageOptimized() for consequences. Definition at line 2655 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::fillComplete | ( | const Teuchos::RCP< Teuchos::ParameterList > & | params = null | ) |
Signal that data entry is complete.
Off-node entries are distributed (via globalAssemble()), repeated entries are summed, and global indices are transformed to local indices.
Definition at line 2628 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::expertStaticFillComplete | ( | const Teuchos::RCP< const map_type > & | domainMap, |
| const Teuchos::RCP< const map_type > & | rangeMap, | ||
| const Teuchos::RCP< const import_type > & | importer = Teuchos::null, |
||
| const Teuchos::RCP< const export_type > & | exporter = Teuchos::null, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = Teuchos::null |
||
| ) |
Perform a fillComplete on a graph that already has data, via setAllIndices().
The graph must already have filled local 1-D storage (lclInds1D_ and rowPtrs_). If the graph has been constructed in any other way, this method will throw an exception. This routine is needed to support other Trilinos packages and should not be called by ordinary users.
Definition at line 2772 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Comm< int > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getComm | ( | ) | const [virtual] |
Returns the communicator.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 483 of file Tpetra_CrsGraph_def.hpp.
| RCP< Node > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNode | ( | ) | const [virtual] |
Returns the underlying node.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 344 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowMap | ( | ) | const [virtual] |
Returns the Map that describes the row distribution in this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 351 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getColMap | ( | ) | const [virtual] |
Returns the Map that describes the column distribution in this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 358 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getDomainMap | ( | ) | const [virtual] |
Returns the Map associated with the domain of this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 365 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRangeMap | ( | ) | const [virtual] |
Returns the Map associated with the domain of this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 372 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getImporter | ( | ) | const [virtual] |
Returns the importer associated with this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 379 of file Tpetra_CrsGraph_def.hpp.
| RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getExporter | ( | ) | const [virtual] |
Returns the exporter associated with this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 386 of file Tpetra_CrsGraph_def.hpp.
| global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumRows | ( | ) | const [virtual] |
Returns the number of global rows in the graph.
Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 297 of file Tpetra_CrsGraph_def.hpp.
| global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumCols | ( | ) | const [virtual] |
Returns the number of global columns in the graph.
Returns the number of entries in the domain map of the matrix. Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 303 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeNumRows | ( | ) | const [virtual] |
Returns the number of graph rows owned on the calling node.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 314 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeNumCols | ( | ) | const [virtual] |
Returns the number of columns connected to the locally owned rows of this graph.
Throws std::runtime_error if hasColMap() == false
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 320 of file Tpetra_CrsGraph_def.hpp.
| GlobalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getIndexBase | ( | ) | const [virtual] |
Returns the index base for global indices for this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 489 of file Tpetra_CrsGraph_def.hpp.
| global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumEntries | ( | ) | const [virtual] |
Returns the global number of entries in the graph.
Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 416 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeNumEntries | ( | ) | const [virtual] |
Returns the local number of entries in the graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 422 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow | ( | GlobalOrdinal | globalRow | ) | const [virtual] |
Returns the current number of entries on this node in the specified global row.
Returns OrdinalTraits<size_t>::invalid() if the specified global row does not belong to this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1622 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow | ( | LocalOrdinal | localRow | ) | const [virtual] |
Returns the current number of entries on this node in the specified local row.
Returns OrdinalTraits<size_t>::invalid() if the specified local row is not valid for this graph.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1638 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeAllocationSize | ( | ) | const |
Returns the total number of indices allocated for the graph, across all rows on this node.
This is the allocation available to the user. Actual allocation may be larger, for example, after calling fillComplete(), and thus this does not necessarily reflect the memory consumption of the this graph.
This quantity is computed during the actual allocation. Therefore, if indicesAreAllocated() == false, this method returns OrdinalTraits<size_t>::invalid().
Definition at line 476 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumAllocatedEntriesInGlobalRow | ( | GlobalOrdinal | globalRow | ) | const |
Returns the current number of allocated entries for this node in the specified global row .
Throws exception std::runtime_error if the specified global row does not belong to this node.
Definition at line 1652 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumAllocatedEntriesInLocalRow | ( | LocalOrdinal | localRow | ) | const |
Returns the current number of allocated entries on this node in the specified local row.
Throws exception std::runtime_error if the specified local row is not valid for this node.
Definition at line 1668 of file Tpetra_CrsGraph_def.hpp.
| global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumDiags | ( | ) | const [virtual] |
Returns the number of global diagonal entries, based on global row/column index comparisons.
Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 337 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeNumDiags | ( | ) | const [virtual] |
Returns the number of local diagonal entries, based on global row/column index comparisons.
Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 331 of file Tpetra_CrsGraph_def.hpp.
| global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalMaxNumRowEntries | ( | ) | const [virtual] |
Maximum number of entries in all rows over all processes.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 428 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeMaxNumRowEntries | ( | ) | const [virtual] |
Maximum number of entries in all rows owned by the calling process.
Undefined if isFillActive().
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 434 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::hasColMap | ( | ) | const [virtual] |
Whether the graph has a column Map.
A CrsGraph has a column Map either because it was given to its constructor, or because it was constructed in fillComplete(). Calling fillComplete() always makes a column Map if the graph does not already have one.
A column Map lets the graph
The latter is mainly useful for a graph associated with a CrsMatrix.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 392 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isLowerTriangular | ( | ) | const [virtual] |
Whether the graph is locally lower triangular.
! isFillActive(). If fill is active, this method's behavior is undefined.Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 458 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isUpperTriangular | ( | ) | const [virtual] |
Whether the graph is locally upper triangular.
! isFillActive(). If fill is active, this method's behavior is undefined.Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 452 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isLocallyIndexed | ( | ) | const [virtual] |
Whether column indices are stored using local indices on the calling process.
If the graph stores column indices as local indices on the calling process, this returns true. Otherwise, it returns false. Exactly one of the following is true:
! isLocallyIndexed() && ! isGloballyIndexed() ! isLocallyIndexed() && isGloballyIndexed() isLocallyIndexed() && ! isGloballyIndexed() The first condition means that the calling process does not own any graph entries. The second condition means that the graph does not (yet) have a column Map, so no process can be locally indexed. The third condition means that the graph has a column Map.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 464 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isGloballyIndexed | ( | ) | const [virtual] |
Whether column indices are stored using global indices on the calling process.
If the graph stores column indices as global indices on the calling process, this returns true. Otherwise, it returns false. Exactly one of the following is true:
! isLocallyIndexed() && ! isGloballyIndexed() ! isLocallyIndexed() && isGloballyIndexed() isLocallyIndexed() && ! isGloballyIndexed() The first condition means that the calling process does not own any graph entries. The second condition means that the graph does not (yet) have a column Map, so no process can be locally indexed. The third condition means that the graph has a column Map.
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 470 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isFillComplete | ( | ) | const [virtual] |
Whether fillComplete() has been called and the graph is in compute mode.
This state is the same on all processes in the graph's communicator. Exactly one of the following is true:
Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 440 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isFillActive | ( | ) | const |
Whether resumeFill() has been called and the graph is in edit mode.
This method returns true under the following conditions:
This state is the same on all processes in the graph's communicator. Exactly one of the following is true:
Definition at line 446 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isSorted | ( | ) | const |
Whether graph indices in all rows are known to be sorted.
A fill-complete graph is always sorted, as is a newly constructed graph. A graph is sorted immediately after calling resumeFill(), but any changes to the graph may result in the sorting status becoming unknown (and therefore, presumed unsorted).
Definition at line 513 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isStorageOptimized | ( | ) | const |
Returns true if storage has been optimized.
Optimized storage means that the allocation of each row is equal to the number of entries. The effect is that a pass through the matrix, i.e., during a mat-vec, requires minimal memory traffic. One limitation of optimized storage is that no new indices can be added to the graph.
Definition at line 398 of file Tpetra_CrsGraph_def.hpp.
| ProfileType Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getProfileType | ( | ) | const |
Returns true if the graph was allocated with static data structures.
Definition at line 410 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy | ( | GlobalOrdinal | GlobalRow, |
| const Teuchos::ArrayView< GlobalOrdinal > & | Indices, | ||
| size_t & | NumIndices | ||
| ) | const [virtual] |
Get a copy of the column indices (as global indices) in the given row.
| GlobalRow | [in] Global index of the row. |
| Indices | [out] Global column indices owned by the calling process in row GlobalRow. This is a view of an array that you, the caller, must allocate in advance. You may call getNumEntriesInGlobalRow() to find the number of entries in this row on the calling process, or getNumAllocatedEntriesInGlobalRow() to get the maximum number of entries in all rows on the calling process. (The latter may be useful when looping over all rows on the calling process.) |
| NumIndices | [out] The number of column indices stored in Indices. If GlobalRow does not belong to (i.e., is not in the row Map on) the calling process, then Indices is unchanged and NumIndices is returned as Teuchos::OrdinalTraits<size_t>::invalid(). |
This method only returns the column indices owned by the calling process. This is a local operation with respect to MPI interprocess communication. It could be the case that multiple processes store different sets of column indices in the same row. In that case, this method only returns the column indices stored on the calling process. Use the row Map to determine which process(es) own which rows.
This method throws std::runtime_error if the output array Indices is not large enough to hold the column indices in the given row.
You may call this method at any time, whether or not the graph has a column Map or is fill complete.
NumIndices will go away. It will be the caller's responsibility to test whether this is less than the size of Indices. Second, instead of throwing an exception if GlobalRow is not owned by the calling process, this method will simply return zero. (This indicates that the calling process does not own any column indices in that row, which is true even if the row index is not in the row Map on any process.) Third, this method will take a Kokkos::View instead of a Teuchos::ArrayView. Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1769 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy | ( | LocalOrdinal | LocalRow, |
| const Teuchos::ArrayView< LocalOrdinal > & | indices, | ||
| size_t & | NumIndices | ||
| ) | const [virtual] |
Get a copy of the column indices (as local indices) in the given row.
| LocalRow | [in] Local index of the row. |
| indices | [out] Local column indices owned by the calling process in row LocalRow. This is a view of an array that you, the caller, must allocate in advance. You may call getNumEntriesInLocalRow() to find the number of entries in this row on the calling process, or getNumAllocatedEntriesInLocalRow() to get the maximum number of entries in all rows on the calling process. (The latter may be useful when looping over all rows on the calling process.) |
| NumIndices | [out] The number of column indices stored in indices. If LocalRow does not belong to (i.e., is not in the row Map on) the calling process, then indices is unchanged and NumIndices is returned as Teuchos::OrdinalTraits<size_t>::invalid(). |
This method only returns the column indices owned by the calling process. This is a local operation with respect to MPI interprocess communication. It could be the case that multiple processes store different sets of column indices in the same row. In that case, this method only returns the column indices stored on the calling process. Use the row Map to determine which process(es) own which rows.
This method throws std::runtime_error if the output array indices is not large enough to hold the column indices in the given row.
isLocallyIndexed() || hasColMap().NumIndices will go away. It will be the caller's responsibility to test whether this is less than the size of indices. Second, instead of throwing an exception if LocalRow is not owned by the calling process, this method will simply return zero. (This indicates that the calling process does not own any column indices in that row, which is true even if the row index is not in the row Map on any process.) Third, this method will take a Kokkos::View instead of a Teuchos::ArrayView. Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 1702 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView | ( | GlobalOrdinal | GlobalRow, |
| Teuchos::ArrayView< const GlobalOrdinal > & | Indices | ||
| ) | const |
Return a const, nonpersisting view of global indices in the given row.
| GlobalRow | [in] Global index of the row. |
| Indices | [out] View of the column indices (as global indices) owned by the calling process in row GlobalRow. If GlobalRow does not belong to the calling process, then Indices is set to null. |
! isLocallyIndexed () indices.size () == getNumEntriesInGlobalRow (GlobalRow) GlobalRow is not owned by the calling process, this method will simply return an empty view. (This indicates that the calling process does not own any column indices in that row, which is true even if the row index is not in the row Map on any process.) Second, this method will return a Kokkos::View instead of a Teuchos::ArrayView. Definition at line 1853 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView | ( | LocalOrdinal | LocalRow, |
| Teuchos::ArrayView< const LocalOrdinal > & | indices | ||
| ) | const |
Return a const, nonpersisting view of local indices in the given row.
| LocalRow | [in] Local index of the row. |
| Indices | [out] View of the column indices (as local indices) owned by the calling process in row LocalRow. If LocalRow does not belong to the calling process, then indices is set to null. |
! isGloballyIndexed () indices.size () == getNumEntriesInLocalRow (LocalRow) LocalRow is not owned by the calling process, this method will simply return an empty view. (This indicates that the calling process does not own any column indices in that row, which is true even if the row index is not in the row Map on any process.) Second, this method will return a Kokkos::View instead of a Teuchos::ArrayView. Definition at line 1816 of file Tpetra_CrsGraph_def.hpp.
| std::string Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::description | ( | ) | const [virtual] |
Return a simple one-line description of this object.
Reimplemented from Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 3899 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::describe | ( | Teuchos::FancyOStream & | out, |
| const Teuchos::EVerbosityLevel | verbLevel = Teuchos::Describable::verbLevel_default |
||
| ) | const [virtual] |
Print the object with some verbosity level to an FancyOStream object.
Reimplemented from Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 3922 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::checkSizes | ( | const SrcDistObject & | source | ) | [virtual] |
Compare the source and target (this) objects for compatibility.
Implements Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 4045 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute | ( | const SrcDistObject & | source, |
| size_t | numSameIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteToLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteFromLIDs | ||
| ) | [virtual] |
Perform copies and permutations that are local to this process.
| source | [in] On entry, the source object, from which we are distributing. We distribute to the destination object, which is *this object. |
| numSameIDs | [in] The umber of elements that are the same on the source and destination (this) objects. These elements are owned by the same process in both the source and destination objects. No permutation occurs. |
| numPermuteIDs | [in] The number of elements that are locally permuted between the source and destination objects. |
| permuteToLIDs | [in] List of the elements that are permuted. They are listed by their LID in the destination object. |
| permuteFromLIDs | [in] List of the elements that are permuted. They are listed by their LID in the source object. |
Implements Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 4059 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare | ( | const SrcDistObject & | source, |
| const Teuchos::ArrayView< const LocalOrdinal > & | exportLIDs, | ||
| Teuchos::Array< GlobalOrdinal > & | exports, | ||
| const Teuchos::ArrayView< size_t > & | numPacketsPerLID, | ||
| size_t & | constantNumPackets, | ||
| Distributor & | distor | ||
| ) | [virtual] |
Perform any packing or preparation required for communication.
| source | [in] Source object for the redistribution. |
| exportLIDs | [in] List of the entries (as local IDs in the source object) we will be sending to other images. |
| exports | [out] On exit, the buffer for data to send. |
| numPacketsPerLID | [out] On exit, the implementation of this method must do one of two things: set numPacketsPerLID[i] to contain the number of packets to be exported for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation need not fill numPacketsPerLID. |
| constantNumPackets | [out] On exit, 0 if numPacketsPerLID has variable contents (different size for each LID). If nonzero, then it is expected that the number of packets per LID is constant, and that constantNumPackets is that value. |
| distor | [in] The Distributor object we are using. |
Implements Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 4159 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::pack | ( | const Teuchos::ArrayView< const LocalOrdinal > & | exportLIDs, |
| Teuchos::Array< GlobalOrdinal > & | exports, | ||
| const Teuchos::ArrayView< size_t > & | numPacketsPerLID, | ||
| size_t & | constantNumPackets, | ||
| Distributor & | distor | ||
| ) | const [virtual] |
Pack this object's data for Import or Export.
Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 4188 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine | ( | const Teuchos::ArrayView< const LocalOrdinal > & | importLIDs, |
| const Teuchos::ArrayView< const GlobalOrdinal > & | imports, | ||
| const Teuchos::ArrayView< size_t > & | numPacketsPerLID, | ||
| size_t | constantNumPackets, | ||
| Distributor & | distor, | ||
| CombineMode | CM | ||
| ) | [virtual] |
Perform any unpacking and combining after communication.
| importLIDs | [in] List of the entries (as LIDs in the destination object) we received from other images. |
| imports | [in] Buffer containing data we received. |
| numPacketsPerLID | [in] If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i]. |
| constantNumPackets | [in] If nonzero, then numPacketsPerLID is constant (same value in all entries) and constantNumPackets is that value. If zero, then numPacketsPerLID[i] is the number of packets imported for importLIDs[i]. |
| distor | [in] The Distributor object we are using. |
| CM | [in] The combine mode to use when combining the imported entries with existing entries. |
Implements Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.
Definition at line 4246 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesPerLocalRowUpperBound | ( | Teuchos::ArrayRCP< const size_t > & | boundPerLocalRow, |
| size_t & | boundForAllLocalRows, | ||
| bool & | boundSameForAllLocalRows | ||
| ) | const |
Get an upper bound on the number of entries that can be stored in each row.
When a CrsGraph is constructed, callers must give an upper bound on the number of entries in each local row. They may either supply a single integer which is the upper bound for all local rows, or they may give an array with a possibly different upper bound for each local row.
This method returns the upper bound for each row. If numEntriesPerLocalRowBound is Teuchos::null on output and boundSameForAllLocalRows is true on output, then numEntriesAllLocalRowsBound is the upper bound for all local rows. If boundSameForAllLocalRows is false on output, then numEntriesPerLocalRowBound has zero or more entries on output, and numEntriesPerLocalRowBound[i_local] is the upper bound for local row i_local.
The output argument boundSameForAllLocalRows is conservative; it only tells us whether boundForAllLocalRows has a meaningful value on output. We don't necessarily check whether all entries of boundPerLocalRow are the same.
Definition at line 2220 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setAllIndices | ( | const Teuchos::ArrayRCP< size_t > & | rowPointers, |
| const Teuchos::ArrayRCP< LocalOrdinal > & | columnIndices | ||
| ) |
Set the graph's data directly, using 1-D storage.
hasColMap() == true rowPointers.size() != getNodeNumRows()+1 Definition at line 2196 of file Tpetra_CrsGraph_def.hpp.
| ArrayRCP< const size_t > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodeRowPtrs | ( | ) | const |
Get an ArrayRCP of the row-offsets.
The returned buffer exists in host-memory. This method may return Teuchos::null if "Delete Row Pointers" was true on fillComplete().
Definition at line 1682 of file Tpetra_CrsGraph_def.hpp.
| ArrayRCP< const LocalOrdinal > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNodePackedIndices | ( | ) | const |
Get an ArrayRCP of the packed column-indices.
The returned buffer exists in host-memory.
Definition at line 1691 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::replaceColMap | ( | const Teuchos::RCP< const map_type > & | newColMap | ) |
Replace the graph's current column Map with the given Map.
This only replaces the column Map. It does not change the graph's current column indices, or otherwise apply a permutation. For example, suppose that before calling this method, the calling process owns a row containing local column indices [0, 2, 4]. These indices do not change, nor does their order change, as a result of calling this method.
| newColMap | [in] New column Map. Must be nonnull. |
Definition at line 2935 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::reindexColumns | ( | const Teuchos::RCP< const map_type > & | newColMap, |
| const Teuchos::RCP< const import_type > & | newImport = Teuchos::null, |
||
| const bool | sortIndicesInEachRow = true |
||
| ) |
Reindex the column indices in place, and replace the column Map. Optionally, replace the Import object as well.
| newColMap | [in] New column Map. Must be nonnull. |
| newImport | [in] New Import object. Optional; computed if not provided or if null. Computing an Import is expensive, so it is worth providing this if you can. |
| sortIndicesInEachRow | [in] If true, sort the indices in each row after reindexing. |
Definition at line 2953 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::replaceDomainMapAndImporter | ( | const Teuchos::RCP< const map_type > & | newDomainMap, |
| const Teuchos::RCP< const import_type > & | newImporter | ||
| ) |
Replace the current domain Map and Import with the given parameters.
isFillComplete() == true isFillActive() == false Either the given Import object is null, or the target Map of the given Import is the same as this graph's column Map. Either the given Import object is null, or the source Map of the given Import is the same as this graph's domain Map. Definition at line 3219 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace | ( | const Teuchos::RCP< const map_type > & | newMap | ) | [virtual] |
Remove processes owning zero rows from the Maps and their communicator.
| newMap | [in] This must be the result of calling the removeEmptyProcesses() method on the row Map. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes. |
This method satisfies the strong exception guarantee, as long the destructors of Export, Import, and Map do not throw exceptions. This means that either the method returns normally (without throwing an exception), or there are no externally visible side effects. However, this does not guarantee no deadlock when the graph's original communicator contains more than one process. In order to prevent deadlock, you must still wrap this call in a try/catch block and do an all-reduce over all processes in the original communicator to test whether the call succeeded. This safety measure should usually be unnecessary, since the method call should only fail on user error or failure to allocate memory.
Definition at line 4296 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::makeColMap | ( | ) | [protected] |
Make the graph's column Map, if it does not already have one.
Definition at line 3508 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertIndices | ( | const RowInfo & | rowInfo, |
| const SLocalGlobalViews & | newInds, | ||
| const ELocalGlobal | lg, | ||
| const ELocalGlobal | I | ||
| ) | [protected] |
Insert indices into the given row.
! (lg == LocalIndices && I == GlobalIndices). It does not make sense to give this method local column indices (meaning that the graph has a column Map), yet to ask it to store global indices.| rowInfo | [in] Result of CrsGraph's getRowInfo() or updateAllocAndValues() methods, for the locally owned row (whose local index is rowInfo.localRow) for which you want to insert indices. |
| newInds | [in] View of the column indices to insert. If lg == GlobalIndices, then newInds.ginds, a Teuchos::ArrayView<const GlobalOrdinal>, contains the (global) column indices to insert. Otherwise, if lg == LocalIndices, then newInds.linds, a Teuchos::ArrayView<const LocalOrdinal>, contains the (local) column indices to insert. |
| lg | If lg == GlobalIndices, then the input indices (in newInds) are global indices. Otherwise, if lg == LocalIndices, the input indices are local indices. |
| I | If lg == GlobalIndices, then this method will store the input indices as global indices. Otherwise, if I == LocalIndices, this method will store the input indices as local indices. |
Definition at line 1037 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertIndicesAndValues | ( | const RowInfo & | rowInfo, |
| const SLocalGlobalViews & | newInds, | ||
| const Teuchos::ArrayView< Scalar > & | oldRowVals, | ||
| const Teuchos::ArrayView< const Scalar > & | newRowVals, | ||
| const ELocalGlobal | lg, | ||
| const ELocalGlobal | I | ||
| ) | [inline, protected] |
Insert indices and their values into the given row.
| Scalar | The type of a single value. When this method is called by CrsMatrix, Scalar corresponds to the first template parameter of CrsMatrix. |
! (lg == LocalIndices && I == GlobalIndices). It does not make sense to give this method local column indices (meaning that the graph has a column Map), yet to ask it to store global indices.| rowInfo | [in] Result of CrsGraph's getRowInfo() or updateAllocAndValues() methods, for the locally owned row (whose local index is rowInfo.localRow) for which you want to insert indices. |
| newInds | [in] View of the column indices to insert. If lg == GlobalIndices, then newInds.ginds, a Teuchos::ArrayView<const GlobalOrdinal>, contains the (global) column indices to insert. Otherwise, if lg == LocalIndices, then newInds.linds, a Teuchos::ArrayView<const LocalOrdinal>, contains the (local) column indices to insert. |
| oldRowVals | [out] View of the current values. They will be overwritten with the new values. |
| newRowVals | [in] View of the new values. They will be copied over the old values. |
| lg | If lg == GlobalIndices, then the input indices (in newInds) are global indices. Otherwise, if lg == LocalIndices, the input indices are local indices. |
| I | If lg == GlobalIndices, then this method will store the input indices as global indices. Otherwise, if I == LocalIndices, this method will store the input indices as local indices. |
Definition at line 1288 of file Tpetra_CrsGraph_decl.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertLocalIndicesFiltered | ( | const LocalOrdinal | localRow, |
| const Teuchos::ArrayView< const LocalOrdinal > & | indices | ||
| ) | [protected] |
Like insertLocalIndices(), but with column Map filtering.
Definition at line 1974 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndicesFiltered | ( | const GlobalOrdinal | localRow, |
| const Teuchos::ArrayView< const GlobalOrdinal > & | indices | ||
| ) | [protected] |
Like insertGlobalIndices(), but with column Map filtering.
Definition at line 2110 of file Tpetra_CrsGraph_def.hpp.
| LocalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::transformLocalValues | ( | RowInfo | rowInfo, |
| const Teuchos::ArrayView< Scalar > & | rowVals, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | inds, | ||
| const Teuchos::ArrayView< const Scalar > & | newVals, | ||
| BinaryFunction | f | ||
| ) | const [inline, protected] |
Transform the given values using local indices.
| rowInfo | [in] Information about a given row of the graph. |
| rowVals | [in/out] The values to be transformed. They correspond to the row indicated by rowInfo. |
| inds | [in] The (local) indices in the row, for which to transform the corresponding values in rowVals. |
| newVals | [in] Values to use for transforming rowVals. It's probably OK for these to alias rowVals. |
| f | [in] A binary function used to transform rowVals. |
This method transforms the values using the expression
newVals[k] = f( rowVals[k], newVals[j] );
where k is the local index corresponding to inds[j]. It ignores elements of inds that are not owned by the calling process.
Definition at line 1347 of file Tpetra_CrsGraph_decl.hpp.
| LocalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::transformGlobalValues | ( | RowInfo | rowInfo, |
| const Teuchos::ArrayView< Scalar > & | rowVals, | ||
| const Teuchos::ArrayView< const GlobalOrdinal > & | inds, | ||
| const Teuchos::ArrayView< const Scalar > & | newVals, | ||
| BinaryFunction | f | ||
| ) | const [inline, protected] |
Transform the given values using global indices.
| rowInfo | [in] Information about a given row of the graph. |
| rowVals | [in/out] The values to be transformed. They correspond to the row indicated by rowInfo. |
| inds | [in] The (global) indices in the row, for which to transform the corresponding values in rowVals. |
| newVals | [in] Values to use for transforming rowVals. It's probably OK for these to alias rowVals. |
| f | [in] A binary function used to transform rowVals. |
Definition at line 1395 of file Tpetra_CrsGraph_decl.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isMerged | ( | ) | const [protected] |
Whether duplicate column indices in each row have been merged.
Definition at line 522 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setLocallyModified | ( | ) | [protected] |
Report that we made a local modification to its structure.
Call this after making a local change to the graph's structure. Changing the structure locally invalidates the "is sorted" and "is merged" states.
Definition at line 530 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortAllIndices | ( | ) | [protected] |
Sort the column indices in all the rows.
Definition at line 3491 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortRowIndices | ( | RowInfo | rowinfo | ) | [protected] |
Sort the column indices in the given row.
Definition at line 1255 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortRowIndicesAndValues | ( | const RowInfo | rowinfo, |
| const Teuchos::ArrayView< Scalar > & | values | ||
| ) | [protected] |
Sort the column indices and their values in the given row.
| Scalar | The type of the values. When calling this method from CrsMatrix, this should be the same as the Scalar template parameter of CrsMatrix. |
| rowinfo | [in] Result of getRowInfo() for the row. |
| values | [in/out] On input: values for the given row. If indices is an array of the column indices in the row, then values and indices should have the same number of entries, and indices[k] should be the column index corresponding to values[k]. On output: the same values, but sorted in the same order as the (now sorted) column indices in the row. |
Definition at line 1268 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::mergeAllIndices | ( | ) | [protected] |
Merge duplicate row indices in all of the rows.
isGloballyIndexed() == false.isMerged() == false. That is, this function would normally only be called after calling sortIndices(). Definition at line 3837 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::mergeRowIndices | ( | RowInfo | rowinfo | ) | [protected] |
Merge duplicate row indices in the given row.
isStorageOptimized() == false Definition at line 1281 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::mergeRowIndicesAndValues | ( | RowInfo | rowinfo, |
| const Teuchos::ArrayView< Scalar > & | rowValues | ||
| ) | [protected] |
Merge duplicate row indices in the given row, along with their corresponding values.
This method is only called by CrsMatrix, for a CrsMatrix whose graph is this CrsGraph instance. It is only called when the matrix owns the graph, not when the matrix was constructed with a const graph.
isStorageOptimized() == false | void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setDomainRangeMaps | ( | const Teuchos::RCP< const map_type > & | domainMap, |
| const Teuchos::RCP< const map_type > & | rangeMap | ||
| ) | [protected] |
Set the domain and range Maps, and invalidate the Import and/or Export objects if necessary.
If the domain Map has changed, invalidate the Import object (if there is one). Likewise, if the range Map has changed, invalidate the Export object (if there is one).
Definition at line 1368 of file Tpetra_CrsGraph_def.hpp.
| RowInfo Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowInfo | ( | size_t | myRow | ) | const [protected] |
Return indexing information for the given (local) row index.
| myRow | [in] Local row index, as a size_t. |
Definition at line 808 of file Tpetra_CrsGraph_def.hpp.
| ArrayView< const LocalOrdinal > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalView | ( | RowInfo | rowinfo | ) | const [protected] |
Get a const view of the local column indices in the given row.
This is one of the methods that take a RowInfo object returned by getRowInfo(). Note that it returns a view of all space for column indices in the row. The valid entries are the first rowinfo.numEntries entries of the returned view.
Definition at line 730 of file Tpetra_CrsGraph_def.hpp.
| ArrayView< LocalOrdinal > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalViewNonConst | ( | RowInfo | rowinfo | ) | [protected] |
Get a nonconst view of the local column indices in the given row.
This is one of the methods that take a RowInfo object returned by getRowInfo(). Note that it returns a view of all space for column indices in the row. The valid entries are the first rowinfo.numEntries entries of the returned view.
Definition at line 750 of file Tpetra_CrsGraph_def.hpp.
| ArrayView< const GlobalOrdinal > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalView | ( | RowInfo | rowinfo | ) | const [protected] |
Get a const view of the global column indices in the given row.
This is one of the methods that take a RowInfo object returned by getRowInfo(). Note that it returns a view of all space for column indices in the row. The valid entries are the first rowinfo.numEntries entries of the returned view.
Definition at line 770 of file Tpetra_CrsGraph_def.hpp.
| ArrayView< GlobalOrdinal > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalViewNonConst | ( | RowInfo | rowinfo | ) | [protected] |
Get a nonconst view of the global column indices in the given row.
This is one of the methods that take a RowInfo object returned by getRowInfo(). Note that it returns a view of all space for column indices in the row. The valid entries are the first rowinfo.numEntries entries of the returned view.
Definition at line 790 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::findLocalIndex | ( | RowInfo | rowinfo, |
| LocalOrdinal | ind, | ||
| size_t | hint = 0 |
||
| ) | const [protected] |
Find the column offset corresponding to the given (local) column index.
The name of this method is a bit misleading. It does not actually find the column index. Instead, it takes a local column index ind, and returns the corresponding offset into the raw array of column indices (whether that be 1-D or 2-D storage).
| rowinfo | [in] Result of getRowInfo() for the given row. |
| ind | [in] (Local) column index for which to find the offset. |
| hint | [in] Hint for where to find ind in the column indices for the given row. If colInds is the ArrayView of the (local) column indices for the given row, and if colInds[hint] == ind, then the hint is correct. The hint is ignored if it is out of range (that is, greater than or equal to the number of entries in the given row). |
The hint optimizes for the case of calling this method several times with the same row (as it would be in transformLocalValues) when several index inputs occur in consecutive sequence. This may occur (for example) when there are multiple degrees of freedom per mesh point, and users are handling the assignment of degrees of freedom to global indices manually (rather than letting BlockMap take care of it). In that case, users might choose to assign the degrees of freedom for a mesh point to consecutive global indices. Epetra implements the hint for this reason.
The hint only costs two comparisons (one to check range, and the other to see if the hint was correct), and it can save searching for the indices (which may take a lot more than two comparisons).
Definition at line 1388 of file Tpetra_CrsGraph_def.hpp.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::findLocalIndex | ( | RowInfo | rowinfo, |
| LocalOrdinal | ind, | ||
| Teuchos::ArrayView< const LocalOrdinal > | colInds, | ||
| size_t | hint = 0 |
||
| ) | const [protected] |
Find the column offset corresponding to the given (local) column index, given a view of the (local) column indices.
The name of this method is a bit misleading. It does not actually find the column index. Instead, it takes a local column index ind, and returns the corresponding offset into the raw array of column indices (whether that be 1-D or 2-D storage).
It is best to use this method if you plan to call it several times for the same row, like in transformLocalValues(). In that case, it amortizes the overhead of calling getLocalView().
| rowinfo | [in] Result of getRowInfo() for the given row. |
| ind | [in] (Local) column index for which to find the offset. |
| colInds | [in] View of all the (local) column indices for the given row. |
| hint | [in] Hint for where to find ind in the column indices for the given row. If colInds is the ArrayView of the (local) column indices for the given row, and if colInds[hint] == ind, then the hint is correct. The hint is ignored if it is out of range (that is, greater than or equal to the number of entries in the given row). |
See the documentation of the three-argument version of this method for an explanation and justification of the hint.
| size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::findGlobalIndex | ( | RowInfo | rowinfo, |
| GlobalOrdinal | ind, | ||
| size_t | hint = 0 |
||
| ) | const [protected] |
Find the column offset corresponding to the given (global) column index.
The name of this method is a bit misleading. It does not actually find the column index. Instead, it takes a global column index ind, and returns the corresponding offset into the raw array of column indices (whether that be 1-D or 2-D storage).
Definition at line 1455 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::checkInternalState | ( | ) | const [protected] |
Check whether the graph's state is valid.
This method is useful for debugging. It gets called often in a debug build.
Definition at line 1510 of file Tpetra_CrsGraph_def.hpp.
| bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::hasRowInfo | ( | ) | const [protected] |
Whether it is valid to call getRowInfo().
FIXME (mfh 21 Oct 2013, 28 Sep 2014) As far as I can tell, this should always return true. Why do we need it? It looks like, historically, the graph (by default) "deleted row info" at fillComplete(). This is probably why many CrsGraph methods check whether hasRowInfo() returns true before doing anything. However, I think that nonintuitive behavior was fixed later, such that hasRowInfo() should always return true. If this is actually the case, then we should dispense with this method.
Definition at line 4395 of file Tpetra_CrsGraph_def.hpp.
| void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doImport | ( | const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Import< LocalOrdinal, GlobalOrdinal, Node > & | importer, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Import data into this object using an Import object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Import object if you want to do an Import, else use doExport() with a precomputed Export object.
| source | [in] The "source" object for redistribution. |
| importer | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap(). |
| CM | [in] How to combine incoming data with the same global index. |
| void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doImport | ( | const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Export< LocalOrdinal, GlobalOrdinal, Node > & | exporter, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Import data into this object using an Export object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doImport() that takes a precomputed Import object in that case.
| source | [in] The "source" object for redistribution. |
| exporter | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.) |
| CM | [in] How to combine incoming data with the same global index. |
| void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doExport | ( | const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Export< LocalOrdinal, GlobalOrdinal, Node > & | exporter, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Export data into this object using an Export object ("forward mode").
The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use this method with your precomputed Export object if you want to do an Export, else use doImport() with a precomputed Import object.
| source | [in] The "source" object for redistribution. |
| exporter | [in] Precomputed data redistribution plan. Its source Map must be the same as the input DistObject's Map, and its target Map must be the same as this->getMap(). |
| CM | [in] How to combine incoming data with the same global index. |
| void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doExport | ( | const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > & | source, |
| const Import< LocalOrdinal, GlobalOrdinal, Node > & | importer, | ||
| CombineMode | CM | ||
| ) | [inherited] |
Export data into this object using an Import object ("reverse mode").
The input DistObject is always the source of the data redistribution operation, and the *this object is always the target.
If you don't know the difference between forward and reverse mode, then you probably want forward mode. Use the version of doExport() that takes a precomputed Export object in that case.
| source | [in] The "source" object for redistribution. |
| importer | [in] Precomputed data redistribution plan. Its target Map must be the same as the input DistObject's Map, and its source Map must be the same as this->getMap(). (Note the difference from forward mode.) |
| CM | [in] How to combine incoming data with the same global index. |
| bool Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::isDistributed | ( | ) | const [inherited] |
Whether this is a globally distributed object.
For a definition of "globally distributed" (and its opposite, "locally replicated"), see the documentation of Map's isDistributed() method.
| virtual Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::getMap | ( | ) | const [inline, virtual, inherited] |
The Map describing the parallel distribution of this object.
Note that some Tpetra objects might be distributed using multiple Map objects. For example, CrsMatrix has both a row Map and a column Map. It is up to the subclass to decide which Map to use when invoking the DistObject constructor.
Definition at line 316 of file Tpetra_DistObject_decl.hpp.
| void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::print | ( | std::ostream & | os | ) | const [inherited] |
Print this object to the given output stream.
We generally assume that all MPI processes can print to the given stream.
| virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace | ( | const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & | newMap | ) | [virtual, inherited] |
Remove processes which contain no elements in this object's Map.
On input, this object is distributed over the Map returned by getMap() (the "original Map," with its communicator, the "original communicator"). The input newMap of this method must be the same as the result of calling getMap()->removeEmptyProcesses(). On processes in the original communicator which contain zero elements ("excluded processes," as opposed to "included processes"), the input newMap must be Teuchos::null (which is what getMap()->removeEmptyProcesses() returns anyway).
On included processes, reassign this object's Map (that would be returned by getMap()) to the input newMap, and do any work that needs to be done to restore correct semantics. On excluded processes, free any data that needs freeing, and do any other work that needs to be done to restore correct semantics.
This method has collective semantics over the original communicator. On exit, the only method of this object which is safe to call on excluded processes is the destructor. This implies that subclasses' destructors must not contain communication operations.
| virtual size_t Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::constantNumberOfPackets | ( | ) | const [protected, virtual, inherited] |
Whether the implementation's instance promises always to have a constant number of packets per LID, and if so, how many packets per LID there are.
If this method returns zero, the instance says that it might possibly have a different number of packets for each LID to send or receive. If it returns nonzero, the instance promises that the number of packets is the same for all LIDs, and that the return value is this number of packets per LID.
The default implementation of this method returns zero. This does not affect the behavior of doTransfer() in any way. If a nondefault implementation returns nonzero, doTransfer() will use this information to avoid unnecessary allocation and / or resizing of arrays.
| virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doTransfer | ( | const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > & | src, |
| CombineMode | CM, | ||
| size_t | numSameIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteToLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | permuteFromLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | remoteLIDs, | ||
| const Teuchos::ArrayView< const LocalOrdinal > & | exportLIDs, | ||
| Distributor & | distor, | ||
| ReverseOption | revOp | ||
| ) | [protected, virtual, inherited] |
Redistribute data across memory images.
| src | [in] The source object, to redistribute into the target object, which is *this object. |
| CM | [in] The combine mode that describes how to combine values that map to the same global ID on the same process. |
| permuteToLIDs | [in] See copyAndPermute(). |
| permuteFromLIDs | [in] See copyAndPermute(). |
| remoteLIDs | [in] List of entries (as local IDs) in the destination object to receive from other processes. |
| exportLIDs | [in] See packAndPrepare(). |
| distor | [in/out] The Distributor object that knows how to redistribute data. |
| revOp | [in] Whether to do a forward or reverse mode redistribution. |
| virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::createViews | ( | ) | const [protected, virtual, inherited] |
Hook for creating a const view.
doTransfer() calls this on the source object. By default, it does nothing, but the source object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.
| virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::createViewsNonConst | ( | KokkosClassic::ReadWriteOption | rwo | ) | [protected, virtual, inherited] |
Hook for creating a nonconst view.
doTransfer() calls this on the destination (*this) object. By default, it does nothing, but the destination object can use this as a hint to fetch data from a compute buffer on an off-CPU device (such as a GPU) into host memory.
| rwo | [in] Whether to create a write-only or a read-and-write view. For Kokkos Node types where compute buffers live in a separate memory space (e.g., in the device memory of a discrete accelerator like a GPU), a write-only view only requires copying from host memory to the compute buffer, whereas a read-and-write view requires copying both ways (once to read, from the compute buffer to host memory, and once to write, back to the compute buffer). |
| virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::releaseViews | ( | ) | const [protected, virtual, inherited] |
Hook for releasing views.
doTransfer() calls this on both the source and destination objects, once it no longer needs to access that object's data. By default, this method does nothing. Implementations may use this as a hint to free host memory which is a view of a compute buffer, once the host memory view is no longer needed. Some implementations may prefer to mirror compute buffers in host memory; for these implementations, releaseViews() may do nothing.
| Teuchos::RCP< CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > createCrsGraph | ( | const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > & | map, |
| size_t | maxNumEntriesPerRow = 0, |
||
| const Teuchos::RCP< Teuchos::ParameterList > & | params = Teuchos::null |
||
| ) | [related] |
Non-member function to create an empty CrsGraph given a row map and a non-zero profile.
Definition at line 1809 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const map_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::rowMap_ [protected] |
The Map describing the distribution of rows of the graph.
Definition at line 1650 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const map_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::colMap_ [protected] |
The Map describing the distribution of columns of the graph.
Definition at line 1652 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const map_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::rangeMap_ [protected] |
The Map describing the range of the (matrix corresponding to the) graph.
Definition at line 1654 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const map_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::domainMap_ [protected] |
The Map describing the domain of the (matrix corresponding to the) graph.
Definition at line 1656 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const import_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::importer_ [protected] |
Teuchos::RCP<const export_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::exporter_ [protected] |
ProfileType Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::pftype_ [protected] |
Whether the graph was allocated with static or dynamic profile.
Definition at line 1683 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<const size_t> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::numAllocPerRow_ [protected] |
The maximum number of entries to allow in each locally owned row, per row.
This is an argument to some of the graph's constructors. Either this or numAllocForAllRows_ is used, but not both.
If this is not set in the constructor, it is allocated temporarily, if necessary, in allocateIndices(). In that same method, it is used to allocate the row offsets array, then discarded (set to null) unconditionally.
Definition at line 1694 of file Tpetra_CrsGraph_decl.hpp.
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::numAllocForAllRows_ [protected] |
The maximum number of entries to allow in each locally owned row.
This is an argument to some of the graph's constructors. Either this or numAllocPerRow_ is used, but not both.
Definition at line 1700 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<LocalOrdinal> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::lclInds1D_ [protected] |
lclInds1D_ are the indices for all rows
Definition at line 1714 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<GlobalOrdinal> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::gblInds1D_ [protected] |
gblInds1D_ are the indices for all rows
Definition at line 1716 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<Teuchos::Array<LocalOrdinal> > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::lclInds2D_ [protected] |
lclInds2D_[r] are the indices for row r.
Definition at line 1732 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<Teuchos::Array<GlobalOrdinal> > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::gblInds2D_ [protected] |
gblInds2D_[r] are the indices for row r.
Definition at line 1735 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::ArrayRCP<size_t> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::numRowEntries_ [protected] |
The number of local entries in each locally owned row.
This is deallocated in fillComplete() if fillComplete()'s "Optimize Storage" parameter is set to true.
Definition at line 1741 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::lowerTriangular_ [protected] |
Whether the graph is locally lower triangular.
Definition at line 1748 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::upperTriangular_ [protected] |
Whether the graph is locally upper triangular.
Definition at line 1750 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::indicesAreSorted_ [protected] |
Whether the graph's indices are sorted in each row, on this process.
Definition at line 1752 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::noRedundancies_ [protected] |
Whether the graph's indices are non-redundant (merged) in each row, on this process.
Definition at line 1755 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::haveLocalConstants_ [protected] |
Whether this process has computed local constants.
Definition at line 1757 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants_ [protected] |
Whether all processes have computed global constants.
Definition at line 1759 of file Tpetra_CrsGraph_decl.hpp.
std::map<GlobalOrdinal, std::vector<GlobalOrdinal> > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::nonlocals_ [protected] |
Nonlocal data given to insertGlobalValues() or sumIntoGlobalValues().
"Nonlocal" means "corresponding to rows not in the row Map on the calling process." These data are communicated and emptied out in globalAssemble(). That method is in turn called, if necessary, in fillComplete().
Definition at line 1768 of file Tpetra_CrsGraph_decl.hpp.
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortGhostsAssociatedWithEachProcessor_ [protected] |
Whether to require makeColMap() (and therefore fillComplete()) to order column Map GIDs associated with each remote process in ascending order.
makeColMap() always groups remote GIDs by process rank, so that all remote GIDs with the same owning rank occur contiguously. By default, it always sorts remote GIDs in increasing order within those groups. This behavior differs from Epetra, which does not sort remote GIDs with the same owning process.
This is true by default, which means "sort remote GIDs." If you don't want to sort (for compatibility with Epetra), call sortGhostColumnGIDsWithinProcessBlock(false).
Definition at line 1797 of file Tpetra_CrsGraph_decl.hpp.
Teuchos::RCP<const Map<LocalOrdinal,GlobalOrdinal,Node> > Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::map_ [protected, inherited] |
The Map over which this object is distributed.
Definition at line 611 of file Tpetra_DistObject_decl.hpp.
Teuchos::Array<GlobalOrdinal > Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::imports_ [protected, inherited] |
Buffer into which packed data are imported (received from other processes).
Definition at line 615 of file Tpetra_DistObject_decl.hpp.
Teuchos::Array<size_t> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::numImportPacketsPerLID_ [protected, inherited] |
Number of packets to receive for each receive operation.
This array is used in Distributor::doPosts() (and doReversePosts()) when starting the ireceive operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Definition at line 627 of file Tpetra_DistObject_decl.hpp.
Teuchos::Array<GlobalOrdinal > Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::exports_ [protected, inherited] |
Buffer from which packed data are exported (sent to other processes).
Definition at line 630 of file Tpetra_DistObject_decl.hpp.
Teuchos::Array<size_t> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::numExportPacketsPerLID_ [protected, inherited] |
Number of packets to send for each send operation.
This array is used in Distributor::doPosts() (and doReversePosts()) for preparing for the send operation.
This may be ignored in doTransfer() if constantNumPackets is nonzero, indicating a constant number of packets per LID. (For example, MultiVector sets the constantNumPackets output argument of packAndPrepare() to the number of columns in the multivector.)
Definition at line 642 of file Tpetra_DistObject_decl.hpp.
1.7.6.1