Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Related Functions | List of all members
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > Class Template Referenceabstract

A distributed graph accessed by rows (adjacency lists) and stored sparsely. More...

#include <Tpetra_CrsGraph_decl.hpp>

Inheritance diagram for Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >:
Inheritance graph
[legend]

Public Types

using local_ordinal_type = LocalOrdinal
 The type of the graph's local indices. More...
 
using global_ordinal_type = GlobalOrdinal
 The type of the graph's global indices. More...
 
using device_type = typename Node::device_type
 This class' Kokkos device type. More...
 
using execution_space = typename device_type::execution_space
 This class' Kokkos execution space. More...
 
using node_type = Node
 This class' Kokkos Node type. More...
 
using local_graph_device_type = KokkosSparse::StaticCrsGraph< local_ordinal_type, Kokkos::LayoutLeft, device_type, void, size_t >
 The type of the part of the sparse graph on each MPI process. More...
 
using local_graph_host_type = typename local_graph_device_type::HostMirror
 The type of the part of the sparse graph on each MPI process. More...
 
using map_type = ::Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node >
 The Map specialization used by this class. More...
 
using import_type = ::Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >
 The Import specialization used by this class. More...
 
using export_type = ::Tpetra::Export< LocalOrdinal, GlobalOrdinal, Node >
 The Export specialization used by this class. More...
 
using local_inds_device_view_type = typename row_graph_type::local_inds_device_view_type
 The Kokkos::View type for views of local ordinals on device and host. More...
 
using global_inds_device_view_type = typename row_graph_type::global_inds_device_view_type
 The Kokkos::View type for views of global ordinals on device and host. More...
 

Public Member Functions

void importAndFillComplete (Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &destGraph, const import_type &importer, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
 Import from this to the given destination graph, and make the result fill complete. More...
 
void importAndFillComplete (Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &destGraph, const import_type &rowImporter, const import_type &domainImporter, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const
 Import from this to the given destination graph, and make the result fill complete. More...
 
void exportAndFillComplete (Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &destGraph, const export_type &exporter, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null) const
 Export from this to the given destination graph, and make the result fill complete. More...
 
void exportAndFillComplete (Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &destGraph, const export_type &rowExporter, const export_type &domainExporter, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params) const
 Export from this to the given destination graph, and make the result fill complete. More...
 
bool haveGlobalConstants () const
 Returns true if globalConstants have been computed; false otherwise. More...
 
void computeGlobalConstants ()
 Compute global constants, if they have not yet been computed. More...
 
local_graph_device_type getLocalGraphDevice () const
 Get the local graph. More...
 
Constructor/Destructor Methods
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a single upper bound for the number of entries in all rows on the calling process. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Kokkos::DualView< const size_t *, device_type > &numEntPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a (possibly different) upper bound for the number of entries in each row. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::ArrayView< const size_t > &numEntPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying a (possibly different) upper bound for the number of entries in each row (legacy KokkosClassic version). More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and a single upper bound for the number of entries in all rows on the calling process. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Kokkos::DualView< const size_t *, device_type > &numEntPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and number of entries in each row. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayView< const size_t > &numEntPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and number of entries in each row (legacy KokkosClassic version). More...
 
 CrsGraph (CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &originalGraph, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and an existing graph to subview. The graph created will point to the views of the existing graph, but only have the rows contained in the passed-in rowMap. This constructor assumes it will alias the first N rows of the graph, where N is the number of rows in rowMap. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const typename local_graph_device_type::row_map_type &rowPointers, const typename local_graph_device_type::entries_type::non_const_type &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and arrays containing the graph. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< local_ordinal_type > &columnIndices, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and arrays containing the graph. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params. More...
 
 CrsGraph (const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const local_graph_device_type &lclGraph, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column Map and a local graph, which the resulting CrsGraph views. In almost all cases the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params. More...
 
 CrsGraph (const local_graph_device_type &lclGraph, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap=Teuchos::null, const Teuchos::RCP< const map_type > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying column, domain and range maps, and a local graph, which the resulting CrsGraph views. In almost all cases the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params. More...
 
 CrsGraph (const local_graph_device_type &lclGraph, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer, const Teuchos::RCP< const export_type > &exporter, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 
 CrsGraph (const row_ptrs_device_view_type &rowPointers, const local_inds_wdv_type &columnIndices, const Teuchos::RCP< const map_type > &rowMap, const Teuchos::RCP< const map_type > &colMap, const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< const import_type > &importer, const Teuchos::RCP< const export_type > &exporter, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Constructor specifying the local row pointer and column indices arrays of the local graph; the row, column, domain, and range maps; and the importer and exporter. In almost all cases, the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params. More...
 
 CrsGraph (const CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &)=default
 Copy constructor (default). More...
 
CrsGraphoperator= (const CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &)=default
 Assignment operator (default). More...
 
 CrsGraph (CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &&)=default
 Move constructor (default). More...
 
CrsGraphoperator= (CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &&)=default
 Move assignment (default). More...
 
virtual ~CrsGraph ()=default
 Destructor (virtual for memory safety of derived classes). More...
 
bool isIdenticalTo (const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > &graph) const
 Create a cloned CrsGraph for a different Node type. More...
 
Implementation of Teuchos::ParameterListAcceptor
void setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &params) override
 Set the given list of parameters (must be nonnull). More...
 
Teuchos::RCP< const
Teuchos::ParameterList > 
getValidParameters () const override
 Default parameter list suitable for validation. More...
 
Insertion/Removal Methods
void insertGlobalIndices (const global_ordinal_type globalRow, const Teuchos::ArrayView< const global_ordinal_type > &indices)
 Insert global indices into the graph. More...
 
void insertGlobalIndices (const global_ordinal_type globalRow, const local_ordinal_type numEnt, const global_ordinal_type inds[])
 Epetra compatibility version of insertGlobalIndices (see above) that takes input as a raw pointer, rather than Teuchos::ArrayView. More...
 
void insertLocalIndices (const local_ordinal_type localRow, const Teuchos::ArrayView< const local_ordinal_type > &indices)
 Insert local indices into the graph. More...
 
void insertLocalIndices (const local_ordinal_type localRow, const local_ordinal_type numEnt, const local_ordinal_type inds[])
 Epetra compatibility version of insertLocalIndices (see above) that takes input as a raw pointer, rather than Teuchos::ArrayView. More...
 
void removeLocalIndices (local_ordinal_type localRow)
 Remove all graph indices from the specified local row. More...
 
Collective methods for changing the graph's global state
void globalAssemble ()
 Communicate nonlocal contributions to other processes. More...
 
void resumeFill (const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Resume fill operations. More...
 
void fillComplete (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Tell the graph that you are done changing its structure. More...
 
void fillComplete (const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Tell the graph that you are done changing its structure; set default domain and range Maps. More...
 
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 > &params=Teuchos::null)
 Perform a fillComplete on a graph that already has data, via setAllIndices(). More...
 
Methods implementing RowGraph.
Teuchos::RCP< const
Teuchos::Comm< int > > 
getComm () const override
 Returns the communicator. More...
 
Teuchos::RCP< const map_typegetRowMap () const override
 Returns the Map that describes the row distribution in this graph. More...
 
Teuchos::RCP< const map_typegetColMap () const override
 Returns the Map that describes the column distribution in this graph. More...
 
Teuchos::RCP< const map_typegetDomainMap () const override
 Returns the Map associated with the domain of this graph. More...
 
Teuchos::RCP< const map_typegetRangeMap () const override
 Returns the Map associated with the domain of this graph. More...
 
Teuchos::RCP< const import_typegetImporter () const override
 Returns the importer associated with this graph. More...
 
Teuchos::RCP< const export_typegetExporter () const override
 Returns the exporter associated with this graph. More...
 
global_size_t getGlobalNumRows () const override
 Returns the number of global rows in the graph. More...
 
global_size_t getGlobalNumCols () const override
 Returns the number of global columns in the graph. More...
 
size_t getLocalNumRows () const override
 Returns the number of graph rows owned on the calling node. More...
 
size_t getLocalNumCols () const override
 Returns the number of columns connected to the locally owned rows of this graph. More...
 
global_ordinal_type getIndexBase () const override
 Returns the index base for global indices for this graph. More...
 
global_size_t getGlobalNumEntries () const override
 Returns the global number of entries in the graph. More...
 
size_t getLocalNumEntries () const override
 The local number of entries in the graph. More...
 
size_t getNumEntriesInGlobalRow (global_ordinal_type globalRow) const override
 Returns the current number of entries on this node in the specified global row. More...
 
size_t getNumEntriesInLocalRow (local_ordinal_type localRow) const override
 Get the number of entries in the given row (local index). More...
 
size_t getLocalAllocationSize () const
 The local number of indices allocated for the graph, over all rows on the calling (MPI) process. More...
 
size_t getNumAllocatedEntriesInGlobalRow (global_ordinal_type globalRow) const
 Current number of allocated entries in the given row on the calling (MPI) process, using a global row index. More...
 
size_t getNumAllocatedEntriesInLocalRow (local_ordinal_type localRow) const
 Current number of allocated entries in the given row on the calling (MPI) process, using a local row index. More...
 
size_t getGlobalMaxNumRowEntries () const override
 Maximum number of entries in any row of the graph, over all processes in the graph's communicator. More...
 
size_t getLocalMaxNumRowEntries () const override
 Maximum number of entries in any row of the graph, on this process. More...
 
bool hasColMap () const override
 Whether the graph has a column Map. More...
 
bool isLocallyIndexed () const override
 Whether the graph's column indices are stored as local indices. More...
 
bool isGloballyIndexed () const override
 Whether the graph's column indices are stored as global indices. More...
 
bool isFillComplete () const override
 Whether fillComplete() has been called and the graph is in compute mode. More...
 
bool isFillActive () const
 Whether resumeFill() has been called and the graph is in edit mode. More...
 
bool isSorted () const
 Whether graph indices in all rows are known to be sorted. More...
 
bool isStorageOptimized () const
 Returns true if storage has been optimized. More...
 
void getGlobalRowCopy (global_ordinal_type gblRow, nonconst_global_inds_host_view_type &gblColInds, size_t &numColInds) const override
 Get a copy of the given row, using global indices. More...
 
void getLocalRowCopy (local_ordinal_type gblRow, nonconst_local_inds_host_view_type &gblColInds, size_t &numColInds) const override
 Get a copy of the given row, using local indices. More...
 
void getGlobalRowView (const global_ordinal_type gblRow, global_inds_host_view_type &gblColInds) const override
 Get a const view of the given global row's global column indices. More...
 
bool supportsRowViews () const override
 Whether this class implements getLocalRowView() and getGlobalRowView() (it does). More...
 
void getLocalRowView (const LocalOrdinal lclRow, local_inds_host_view_type &lclColInds) const override
 Get a const view of the given local row's local column indices. More...
 
Overridden from Teuchos::Describable
std::string description () const override
 Return a one-line human-readable description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
 Print this object to the given output stream with the given verbosity level. More...
 
Advanced methods, at increased risk of deprecation.
void getLocalDiagOffsets (const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
 Get offsets of the diagonal entries in the graph. More...
 
void getLocalOffRankOffsets (offset_device_view_type &offsets) const
 Get offsets of the off-rank entries in the graph. More...
 
void getLocalDiagOffsets (Teuchos::ArrayRCP< size_t > &offsets) const
 Backwards compatibility overload of the above method. More...
 
void setAllIndices (const typename local_graph_device_type::row_map_type &rowPointers, const typename local_graph_device_type::entries_type::non_const_type &columnIndices)
 Set the graph's data directly, using 1-D storage. More...
 
void setAllIndices (const Teuchos::ArrayRCP< size_t > &rowPointers, const Teuchos::ArrayRCP< local_ordinal_type > &columnIndices)
 Set the graph's data directly, using 1-D storage. More...
 
row_ptrs_host_view_type getLocalRowPtrsHost () const
 Get a host view of the packed row offsets. More...
 
row_ptrs_device_view_type getLocalRowPtrsDevice () const
 Get a device view of the packed row offsets. More...
 
local_inds_host_view_type getLocalIndicesHost () const
 Get a host view of the packed column indicies. More...
 
local_inds_device_view_type getLocalIndicesDevice () const
 Get a device view of the packed column indicies. More...
 
void replaceColMap (const Teuchos::RCP< const map_type > &newColMap)
 Replace the graph's current column Map with the given Map. More...
 
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. More...
 
void replaceDomainMap (const Teuchos::RCP< const map_type > &newDomainMap)
 Replace the current domain Map with the given objects. More...
 
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. More...
 
void replaceRangeMap (const Teuchos::RCP< const map_type > &newRangeMap)
 Replace the current Range Map with the given objects. More...
 
void replaceRangeMapAndExporter (const Teuchos::RCP< const map_type > &newRangeMap, const Teuchos::RCP< const export_type > &newExporter)
 Replace the current Range Map and Export with the given parameters. More...
 
virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap) override
 Remove processes owning zero rows from the Maps and their communicator. More...
 
Access to entries in a row
virtual void getGlobalRowCopy (const GlobalOrdinal gblRow, nonconst_global_inds_host_view_type &gblColInds, size_t &numColInds) const =0
 Get a copy of the global column indices in a given row of the graph. More...
 
virtual void getLocalRowCopy (const LocalOrdinal lclRow, nonconst_local_inds_host_view_type &lclColInds, size_t &numColInds) const =0
 Get a copy of the local column indices in a given row of the graph. More...
 
virtual void getLocalRowView (const LocalOrdinal lclRow, local_inds_host_view_type &lclColInds) const =0
 Get a constant, nonpersisting, locally indexed view of the given row of the graph. More...
 
virtual void getGlobalRowView (const GlobalOrdinal gblRow, global_inds_host_view_type &gblColInds) const =0
 Get a const, non-persisting view of the given global row's global column indices, as a Teuchos::ArrayView. More...
 
Public methods for redistributing data
void doImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 Import data into this object using an Import object ("forward mode"). More...
 
void doImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 Import data into this object using an Export object ("reverse mode"). More...
 
void doExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 Export data into this object using an Export object ("forward mode"). More...
 
void doExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 Export data into this object using an Import object ("reverse mode"). More...
 
void beginImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void beginImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const CombineMode CM, const bool restrictedMode=false)
 
bool transferArrived () const
 Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step. More...
 
Attribute accessor methods
bool isDistributed () const
 Whether this is a globally distributed object. More...
 
virtual Teuchos::RCP< const
map_type
getMap () const
 The Map describing the parallel distribution of this object. More...
 
I/O methods
void print (std::ostream &os) const
 Print this object to the given output stream. More...
 
Methods for use only by experts
virtual void removeEmptyProcessesInPlace (const Teuchos::RCP< const map_type > &newMap)
 Remove processes which contain no entries in this object's Map. More...
 
const Details::DistributorActor & getActor ()
 Getter for the DistributorActor. More...
 

Protected Member Functions

void setDomainRangeMaps (const Teuchos::RCP< const map_type > &domainMap, const Teuchos::RCP< const map_type > &rangeMap)
 
void computeLocalConstants ()
 Compute local constants, if they have not yet been computed. More...
 
RowInfo getRowInfo (const local_ordinal_type myRow) const
 Get information about the locally owned row with local index myRow. More...
 
RowInfo getRowInfoFromGlobalRowIndex (const global_ordinal_type gblRow) const
 Get information about the locally owned row with global index gblRow. More...
 
void checkInternalState () const
 Throw an exception if the internal state is not consistent. More...
 
void swap (CrsGraph< local_ordinal_type, global_ordinal_type, Node > &graph)
 Swaps the data from *this with the data and maps from graph. More...
 
const row_ptrs_device_view_type & getRowPtrsUnpackedDevice () const
 Get the unpacked row pointers on device. More...
 
const row_ptrs_host_view_type & getRowPtrsUnpackedHost () const
 Get the unpacked row pointers on host. Lazily make a copy from device. More...
 
const row_ptrs_device_view_type & getRowPtrsPackedDevice () const
 Get the packed row pointers on device. More...
 
const row_ptrs_host_view_type & getRowPtrsPackedHost () const
 Get the packed row pointers on host. Lazily make a copy from device. More...
 
local_inds_dualv_type::t_host::const_type getLocalIndsViewHost (const RowInfo &rowinfo) const
 Get a const, locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow). More...
 
local_inds_dualv_type::t_dev::const_type getLocalIndsViewDevice (const RowInfo &rowinfo) const
 Get a const, locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow). More...
 
global_inds_dualv_type::t_host::const_type getGlobalIndsViewHost (const RowInfo &rowinfo) const
 Get a const, globally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow). More...
 
global_inds_dualv_type::t_dev::const_type getGlobalIndsViewDevice (const RowInfo &rowinfo) const
 Get a const, globally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow). More...
 
local_inds_dualv_type::t_host getLocalIndsViewHostNonConst (const RowInfo &rowinfo)
 Get a ReadWrite locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow). More...
 
virtual size_t constantNumberOfPackets () const
 Whether the implementation's instance promises always to have a constant number of packets per LID (local index), and if so, how many packets per LID there are. More...
 
virtual void doTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
 Redistribute data across (MPI) processes. More...
 
virtual bool reallocArraysForNumPacketsPerLid (const size_t numExportLIDs, const size_t numImportLIDs)
 Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary. More...
 
void beginTransfer (const SrcDistObject &src, const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &transfer, const char modeString[], const ReverseOption revOp, const CombineMode CM, const bool restrictedMode)
 Implementation detail of doTransfer. More...
 
Methods governing changes between global and local indices
void makeColMap (Teuchos::Array< int > &remotePIDs)
 Make and set the graph's column Map. More...
 
std::pair< size_t, std::string > makeIndicesLocal (const bool verbose=false)
 Convert column indices from global to local. More...
 
void makeImportExport (Teuchos::Array< int > &remotePIDs, const bool useRemotePIDs)
 Make the Import and Export objects, if needed. More...
 

Protected Attributes

Teuchos::RCP< const map_typerowMap_
 The Map describing the distribution of rows of the graph. More...
 
Teuchos::RCP< const map_typecolMap_
 The Map describing the distribution of columns of the graph. More...
 
Teuchos::RCP< const map_typerangeMap_
 The Map describing the range of the (matrix corresponding to the) graph. More...
 
Teuchos::RCP< const map_typedomainMap_
 The Map describing the domain of the (matrix corresponding to the) graph. More...
 
Teuchos::RCP< const import_typeimporter_
 The Import from the domain Map to the column Map. More...
 
Teuchos::RCP< const export_typeexporter_
 The Export from the row Map to the range Map. More...
 
size_t nodeMaxNumRowEntries_
 Local maximum of the number of entries in each row. More...
 
global_size_t globalNumEntries_
 Global number of entries in the graph. More...
 
global_size_t globalMaxNumRowEntries_
 Global maximum of the number of entries in each row. More...
 
local_inds_wdv_type lclIndsUnpacked_wdv
 Local ordinals of column indices for all rows Valid when isLocallyIndexed is true If OptimizedStorage, storage is PACKED after fillComplete If not OptimizedStorate, storage is UNPACKED after fillComplete; that is, the views have storage equal to sizes provided in CrsGraph constructor. More...
 
local_inds_wdv_type lclIndsPacked_wdv
 Local ordinals of column indices for all rows Valid when isLocallyIndexed is true Built during fillComplete or non-fillComplete constructors Storage is PACKED after fillComplete that is, the views have storage equal to sizes provided in CrsGraph constructor. More...
 
global_inds_wdv_type gblInds_wdv
 Global ordinals of column indices for all rows. More...
 
Kokkos::View< const size_t
*, device_type >
::host_mirror_type 
k_numAllocPerRow_
 The maximum number of entries to allow in each locally owned row, per row. More...
 
size_t numAllocForAllRows_
 The maximum number of entries to allow in each locally owned row. More...
 
Details::EStorageStatus storageStatus_
 Status of the graph's storage, when not in a fill-complete state. More...
 
bool indicesAreSorted_ = true
 Whether the graph's indices are sorted in each row, on this process. More...
 
bool noRedundancies_ = true
 Whether the graph's indices are non-redundant (merged) in each row, on this process. More...
 
bool haveLocalConstants_ = false
 Whether this process has computed local constants. More...
 
bool haveGlobalConstants_ = false
 Whether all processes have computed global constants. More...
 
nonlocals_type nonlocals_
 Nonlocal data given to insertGlobalIndices. More...
 
bool sortGhostsAssociatedWithEachProcessor_ = true
 Whether to require makeColMap() (and therefore fillComplete()) to order column Map GIDs associated with each remote process in ascending order. More...
 

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 > &params=Teuchos::null)
 Nonmember function to create an empty CrsGraph given a row Map and the max number of entries allowed locally per row. More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &importer, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsGraph constructor that fuses Import and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &rowImporter, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &domainImporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &domainMap, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsGraph constructor that fuses Import and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &exporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsGraph constructor that fuses Export and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &rowExporter, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &domainExporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &domainMap, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsGraph constructor that fuses Export and fillComplete(). More...
 

Implementation of DistObject

using buffer_device_type = typename dist_object_type::buffer_device_type
 Kokkos::Device specialization for communication buffers. More...
 
using packet_type = global_ordinal_type
 Type of each entry of the DistObject communication buffer. More...
 
using padding_type = Details::CrsPadding< local_ordinal_type, global_ordinal_type >
 
virtual bool checkSizes (const SrcDistObject &source) override
 Compare the source and target (this) objects for compatibility. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
 
void applyCrsPadding (const padding_type &padding, const bool verbose)
 
std::unique_ptr< padding_typecomputeCrsPadding (const RowGraph< local_ordinal_type, global_ordinal_type, node_type > &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const bool verbose) const
 
std::unique_ptr< padding_typecomputeCrsPaddingForImports (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const bool verbose) const
 
std::unique_ptr< padding_typecomputePaddingForCrsMatrixUnpack (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< char *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const bool verbose) const
 
void computeCrsPaddingForSameIDs (padding_type &padding, const RowGraph< local_ordinal_type, global_ordinal_type, node_type > &source, const local_ordinal_type numSameIDs) const
 
void computeCrsPaddingForPermutedIDs (padding_type &padding, const RowGraph< local_ordinal_type, global_ordinal_type, node_type > &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs) const
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) override
 
virtual void pack (const Teuchos::ArrayView< const local_ordinal_type > &exportLIDs, Teuchos::Array< global_ordinal_type > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets) const override
 
void packFillActive (const Teuchos::ArrayView< const local_ordinal_type > &exportLIDs, Teuchos::Array< global_ordinal_type > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets) const
 
void packFillActiveNew (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets) const
 
virtual void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode) override
 

Methods for inserting indices or transforming values

static const bool useAtomicUpdatesByDefault
 Whether transformLocalValues should use atomic updates by default. More...
 
size_t insertIndices (RowInfo &rowInfo, const SLocalGlobalViews &newInds, const ELocalGlobal lg, const ELocalGlobal I)
 Insert indices into the given row. More...
 
size_t insertGlobalIndicesImpl (const local_ordinal_type lclRow, const global_ordinal_type inputGblColInds[], const size_t numInputInds)
 Insert global indices, using an input local row index. More...
 
size_t insertGlobalIndicesImpl (const RowInfo &rowInfo, const global_ordinal_type inputGblColInds[], const size_t numInputInds, std::function< void(const size_t, const size_t, const size_t)> fun=std::function< void(const size_t, const size_t, const size_t)>())
 Insert global indices, using an input RowInfo. More...
 
void insertLocalIndicesImpl (const local_ordinal_type lclRow, const Teuchos::ArrayView< const local_ordinal_type > &gblColInds, std::function< void(const size_t, const size_t, const size_t)> fun=std::function< void(const size_t, const size_t, const size_t)>())
 
size_t findGlobalIndices (const RowInfo &rowInfo, const Teuchos::ArrayView< const global_ordinal_type > &indices, std::function< void(const size_t, const size_t, const size_t)> fun) const
 Finds indices in the given row. More...
 
void insertGlobalIndicesFiltered (const local_ordinal_type lclRow, const global_ordinal_type gblColInds[], const local_ordinal_type numGblColInds)
 Like insertGlobalIndices(), but with column Map filtering. More...
 
void insertGlobalIndicesIntoNonownedRows (const global_ordinal_type gblRow, const global_ordinal_type gblColInds[], const local_ordinal_type numGblColInds)
 Implementation of insertGlobalIndices for nonowned rows. More...
 

Methods for sorting and merging column indices.

bool isMerged () const
 Whether duplicate column indices in each row have been merged. More...
 
void setLocallyModified ()
 Report that we made a local modification to its structure. More...
 
size_t sortAndMergeRowIndices (const RowInfo &rowInfo, const bool sorted, const bool merged)
 Sort and merge duplicate column indices in the given row. More...
 

Graph data structures (packed and unpacked storage).

typedef Kokkos::View< size_t
*, Kokkos::LayoutLeft,
device_type >
::host_mirror_type 
num_row_entries_type
 Row offsets for "1-D" storage. More...
 
num_row_entries_type k_numRowEntries_
 The number of local entries in each locally owned row. More...
 
offset_device_view_type k_offRankOffsets_
 The offsets for off-rank entries. More...
 

Methods implemented by subclasses and used by doTransfer().

The doTransfer() method uses the subclass' implementations of these methods to implement data transfer. Subclasses of DistObject must implement these methods. This is an instance of the Template Method Pattern. ("Template" here doesn't mean "C++ template"; it means "pattern with holes that are filled in by the subclass' method implementations.")

Teuchos::RCP< const map_typemap_
 The Map over which this object is distributed. More...
 
Kokkos::DualView< packet_type
*, buffer_device_type
imports_
 Buffer into which packed data are imported (received from other processes). More...
 
Kokkos::DualView< size_t
*, buffer_device_type
numImportPacketsPerLID_
 Number of packets to receive for each receive operation. More...
 
Kokkos::DualView< packet_type
*, buffer_device_type
exports_
 Buffer from which packed data are exported (sent to other processes). More...
 
Kokkos::DualView< size_t
*, buffer_device_type
numExportPacketsPerLID_
 Number of packets to send for each send operation. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM)
 Perform copies and permutations that are local to the calling (MPI) process. More...
 
virtual void copyAndPermute (const SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM, const execution_space &space)
 Same as copyAndPermute, but do operations in space. More...
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
 Pack data and metadata for communication (sends). More...
 
virtual void packAndPrepare (const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets, const execution_space &space)
 Same as packAndPrepare, but in an execution space instance. More...
 
virtual void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
 Perform any unpacking and combining after communication. More...
 
virtual void unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode, const execution_space &space)
 
std::unique_ptr< std::string > createPrefix (const char className[], const char methodName[]) const
 
virtual bool reallocImportsIfNeeded (const size_t newSize, const bool verbose, const std::string *prefix, const bool remoteLIDsContiguous=false, const CombineMode CM=INSERT)
 Reallocate imports_ if needed. More...
 

Detailed Description

template<class LocalOrdinal, class GlobalOrdinal, class Node>
class Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >

A distributed graph accessed by rows (adjacency lists) and stored sparsely.

Template Parameters
LocalOrdinalThe type of local indices. See the documentation of Map for requirements.
GlobalOrdinalThe type of global indices. See the documentation of Map for requirements.
NodeThe 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.

Prerequisites

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.

Local vs. global indices and nonlocal insertion

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 190 of file Tpetra_CrsGraph_decl.hpp.

Member Typedef Documentation

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type = LocalOrdinal

The type of the graph's local indices.

Definition at line 208 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type = GlobalOrdinal

The type of the graph's global indices.

Definition at line 210 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::device_type = typename Node::device_type

This class' Kokkos device type.

Definition at line 212 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::execution_space = typename device_type::execution_space

This class' Kokkos execution space.

Definition at line 214 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::node_type = Node

This class' Kokkos Node type.

This is a leftover that will be deprecated and removed. See e.g., GitHub Issue #57.

Definition at line 220 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_device_type = KokkosSparse::StaticCrsGraph<local_ordinal_type, Kokkos::LayoutLeft, device_type, void, size_t>

The type of the part of the sparse graph on each MPI process.

Definition at line 225 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_host_type = typename local_graph_device_type::HostMirror

The type of the part of the sparse graph on each MPI process.

Definition at line 231 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type = ::Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>

The Map specialization used by this class.

Definition at line 235 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::import_type = ::Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>

The Import specialization used by this class.

Definition at line 237 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::export_type = ::Tpetra::Export<LocalOrdinal, GlobalOrdinal, Node>

The Export specialization used by this class.

Definition at line 239 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_device_view_type = typename row_graph_type::local_inds_device_view_type

The Kokkos::View type for views of local ordinals on device and host.

Definition at line 263 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::global_inds_device_view_type = typename row_graph_type::global_inds_device_view_type

The Kokkos::View type for views of global ordinals on device and host.

Definition at line 271 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::buffer_device_type

Kokkos::Device specialization for communication buffers.

See #1088 for why this is not just device_type::device_type.

Definition at line 1169 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
using Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::packet_type = global_ordinal_type

Type of each entry of the DistObject communication buffer.

Definition at line 1171 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
typedef Kokkos::View<size_t*, Kokkos::LayoutLeft, device_type>::host_mirror_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::num_row_entries_type
protected

Row offsets for "1-D" storage.

This is only allocated if "1-D" storage is active. In that case, if beg = k_rowPtrs_(i_lcl) and end = k_rowPtrs_(i_lcl+1) for local row index i_lcl, then

  • if the graph is locally indexed, k_lclInds1D_(beg:end-1) (inclusive range) is the space for any local column indices in local row i_lcl, else
  • if the graph is globally indexed, k_gblInds1D_(beg:end-1) (inclusive range) is the space for any global column indices in local row i_lcl.

Only the first k_numRowEntries_(i_lcl) of these entries are actual valid column indices. Any remaining entries are "extra space." If the graph's storage is packed, then there is no extra space, and the k_numRowEntries_ array is invalid.

If it is allocated, k_rowPtrs_ has length getLocalNumRows()+1. The k_numRowEntries_ array has has length getLocalNumRows(), again if it is allocated. The type of k_numRowEntries_ (see below).

This View gets used only on host. However, making this literally a host View (of Kokkos::HostSpace) causes inexplicable test failures only on CUDA. Thus, I left it as a host_mirror_type, which means (given Trilinos' current UVM requirement) that it will be a UVM allocation.

Definition at line 2402 of file Tpetra_CrsGraph_decl.hpp.

Constructor & Destructor Documentation

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const size_t  maxNumEntriesPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying a single upper bound for the number of entries in all rows on the calling process.

Parameters
rowMap[in] Distribution of rows of the graph.
maxNumEntriesPerRow[in] Maximum number of graph entries per row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 249 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Kokkos::DualView< const size_t *, device_type > &  numEntPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying a (possibly different) upper bound for the number of entries in each row.

Parameters
rowMap[in] Distribution of rows of the graph.
numEntPerRow[in] Maximum number of graph entries to allocate for each row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 342 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::ArrayView< const size_t > &  numEntPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying a (possibly different) upper bound for the number of entries in each row (legacy KokkosClassic version).

Parameters
rowMap[in] Distribution of rows of the graph.
numEntPerRow[in] Maximum number of graph entries to allocate for each row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 292 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const size_t  maxNumEntriesPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and a single upper bound for the number of entries in all rows on the calling process.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
maxNumEntriesPerRow[in] Maximum number of graph entries per row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 269 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Kokkos::DualView< const size_t *, device_type > &  numEntPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and number of entries in each row.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
numEntPerRow[in] Maximum number of graph entries to allocate for each row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 375 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::ArrayView< const size_t > &  numEntPerRow,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and number of entries in each row (legacy KokkosClassic version).

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
numEntPerRow[in] Maximum number of graph entries to allocate for each row. This is a strict upper bound.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 410 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &  originalGraph,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and an existing graph to subview. The graph created will point to the views of the existing graph, but only have the rows contained in the passed-in rowMap. This constructor assumes it will alias the first N rows of the graph, where N is the number of rows in rowMap.

Parameters
rowMap[in] Distribution of rows of the graph.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 462 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const typename local_graph_device_type::row_map_type &  rowPointers,
const typename local_graph_device_type::entries_type::non_const_type &  columnIndices,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and arrays containing the graph. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
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. Entries in each row must be sorted (by local index).
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 494 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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< local_ordinal_type > &  columnIndices,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and arrays containing the graph. In almost all cases the indices must be sorted on input, but if they aren't sorted, "sorted" must be set to false in params.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
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. Entries in each row must be sorted (by local index).
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 519 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const local_graph_device_type lclGraph,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column Map and a local graph, which the resulting CrsGraph views. In almost all cases the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Unlike most other CrsGraph constructors, successful completion of this constructor will result in a fill-complete graph.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph. See replaceColMap() for the requirements.
lclGraph[in] A locally indexed Kokkos::StaticCrsGraph whose local row indices come from the specified row Map, and whose local column indices come from the specified column Map.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 544 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const local_graph_device_type lclGraph,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap = Teuchos::null,
const Teuchos::RCP< const map_type > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying column, domain and range maps, and a local graph, which the resulting CrsGraph views. In almost all cases the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Unlike most other CrsGraph constructors, successful completion of this constructor will result in a fill-complete graph.

Parameters
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph.
domainMap[in] The graph's domain Map. MUST be one to one!
rangeMap[in] The graph's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
lclGraph[in] A locally indexed Kokkos::StaticCrsGraph whose local row indices come from the specified row Map, and whose local column indices come from the specified column Map.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 557 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const local_graph_device_type lclGraph,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< const import_type > &  importer,
const Teuchos::RCP< const export_type > &  exporter,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Create a fill-complete CrsGraph from all the things it needs.

Parameters
lclGraph[in] The local graph. In almost all cases the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Definition at line 626 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const row_ptrs_device_view_type &  rowPointers,
const local_inds_wdv_type columnIndices,
const Teuchos::RCP< const map_type > &  rowMap,
const Teuchos::RCP< const map_type > &  colMap,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< const import_type > &  importer,
const Teuchos::RCP< const export_type > &  exporter,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Constructor specifying the local row pointer and column indices arrays of the local graph; the row, column, domain, and range maps; and the importer and exporter. In almost all cases, the local graph must be sorted on input, but if it isn't sorted, "sorted" must be set to false in params.

Unlike most other CrsGraph constructors, successful completion of this constructor will result in a fill-complete graph.

Parameters
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. Entries in each row must be sorted (by local index).
rowMap[in] Distribution of rows of the graph.
colMap[in] Distribution of columns of the graph.
domainMap[in] The graph's domain Map. MUST be one to one!
rangeMap[in] The graph's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
importer[in] Import from the graph's domain Map to its column Map. If no Import is necessary (i.e., if the domain and column Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
exporter[in] Export from the graph's row Map to its range Map. If no Export is necessary (i.e., if the row and range Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 678 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( const CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &  )
default

Copy constructor (default).

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::CrsGraph ( CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &&  )
default

Move constructor (default).

template<class LocalOrdinal, class GlobalOrdinal, class Node>
virtual Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::~CrsGraph ( )
virtualdefault

Destructor (virtual for memory safety of derived classes).

Note
To Tpetra developers: See the C++ Core Guidelines C.21 ("If you define or <tt>=delete</tt> any default operation, define or <tt>=delete</tt> them all"), in particular the AbstractBase example, for why this destructor declaration implies that we need the above four =default declarations for copy construction, move construction, copy assignment, and move assignment.

Member Function Documentation

template<class LocalOrdinal, class GlobalOrdinal, class Node>
CrsGraph& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::operator= ( const CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &  )
default

Assignment operator (default).

template<class LocalOrdinal, class GlobalOrdinal, class Node>
CrsGraph& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::operator= ( CrsGraph< local_ordinal_type, global_ordinal_type, node_type > &&  )
default

Move assignment (default).

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isIdenticalTo ( const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > &  graph) const

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.

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

  • "Static profile clone" [bool, default: true] If true, creates the clone with a static allocation profile. If false, a dynamic allocation profile is used.
  • "Locally indexed clone" [bool] If 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.
  • "fillComplete clone" [boolean, default: true] If 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. True if and only if CrsGraph is identical to this CrsGraph
Warning
THIS METHOD IS FOR TPETRA DEVELOPERS ONLY. DO NOT RELY ON THIS METHOD. WE MAKE NO PROMISES OF BACKWARDS COMPATIBILITY.

This performs exact matches on objects with in the graphs. That is, internal data structures such as arrays must match exactly in both content and order. This is not performing any sort of isomorphic search.

Parameters
graph[in] a CrsGraph to compare against this one.
Returns
True if the other CrsGraph's data structure is identical to this CrsGraph.

Definition at line 6960 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  params)
override

Set the given list of parameters (must be nonnull).

Definition at line 766 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Teuchos::ParameterList > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getValidParameters ( ) const
override

Default parameter list suitable for validation.

Definition at line 732 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndices ( const global_ordinal_type  globalRow,
const Teuchos::ArrayView< const global_ordinal_type > &  indices 
)

Insert global indices into the graph.

Precondition
globalRow is a valid index in the row Map. It need not be owned by the calling process.
isLocallyIndexed() == false
isStorageOptimized() == false
Postcondition
indicesAreAllocated() == true
isGloballyIndexed() == true

If 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 2381 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndices ( const global_ordinal_type  globalRow,
const local_ordinal_type  numEnt,
const global_ordinal_type  inds[] 
)

Epetra compatibility version of insertGlobalIndices (see above) that takes input as a raw pointer, rather than Teuchos::ArrayView.

Arguments are the same and in the same order as Epetra_CrsGraph::InsertGlobalIndices.

Definition at line 2309 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertLocalIndices ( const local_ordinal_type  localRow,
const Teuchos::ArrayView< const local_ordinal_type > &  indices 
)

Insert local indices into the graph.

Precondition
localRow is a local row belonging to the graph on this process.
isGloballyIndexed() == false
isStorageOptimized() == false
hasColMap() == true
Postcondition
indicesAreAllocated() == true
isLocallyIndexed() == true
Note
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 at insertion or during the next call to fillComplete().

Definition at line 2234 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertLocalIndices ( const local_ordinal_type  localRow,
const local_ordinal_type  numEnt,
const local_ordinal_type  inds[] 
)

Epetra compatibility version of insertLocalIndices (see above) that takes input as a raw pointer, rather than Teuchos::ArrayView.

Arguments are the same and in the same order as Epetra_CrsGraph::InsertMyIndices.

Definition at line 2300 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::removeLocalIndices ( local_ordinal_type  localRow)

Remove all graph indices from the specified local row.

Precondition
localRow is a local row of this graph.
isGloballyIndexed() == false
isStorageOptimized() == false
Postcondition
getNumEntriesInLocalRow(localRow) == 0
indicesAreAllocated() == true
isLocallyIndexed() == true

Definition at line 2465 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::globalAssemble ( )

Communicate nonlocal contributions to other processes.

This method is called automatically by fillComplete(). Most users do not need to call this themselves.

This method must be called collectively (that is, like any MPI collective) over all processes in the graph's communicator.

Definition at line 2657 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::resumeFill ( const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null)

Resume fill operations.

After calling fillComplete(), resumeFill() must be called before initiating any changes to the graph.

resumeFill() may be called repeatedly.

Warning
A CrsGraph instance does not currently (as of 23 Jul 2017) and never did support arbitrary structure changes after the first fillComplete call on that instance. The safest thing to do is not to change structure at all after first fillComplete.
Postcondition
isFillActive() == true
isFillComplete() == false

This method must be called collectively (that is, like any MPI collective) over all processes in the graph's communicator.

Definition at line 2870 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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 = Teuchos::null 
)

Tell the graph that you are done changing its structure.

This tells the graph to optimize its data structures for computational kernels, and to prepare (MPI) communication patterns.

Off-process indices are distributed (via globalAssemble()), indices are sorted, redundant indices are eliminated, and global indices are transformed to local indices.

This method must be called collectively (that is, like any MPI collective) over all processes in the graph's communicator.

Warning
The domain Map and row Map arguments to this method MUST be one to one! If you have Maps that are not one to one, and you do not know how to make a Map that covers the same global indices but is one to one, then you may call Tpetra::createOneToOne() (see Map's header file) to make a one-to-one version of your Map.
Precondition
isFillActive() && ! isFillComplete()
Postcondition
! isFillActive() && isFillComplete()
Parameters
domainMap[in] The graph's domain Map. MUST be one to one!
rangeMap[in] The graph's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
params[in/out] List of parameters controlling this method's behavior. See below for valid parameters.

List of valid parameters in params:

Definition at line 2905 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::fillComplete ( const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null)

Tell the graph that you are done changing its structure; set default domain and range Maps.

See above three-argument version of fillComplete for full documentation. If the graph does not yet have domain and range Maps (i.e., if fillComplete has not yet been called on this graph at least once), then this method uses the graph's row Map (result of this->getRowMap()) as both the domain Map and the range Map. Otherwise, this method uses the graph's existing domain and range Maps.

This method must be called collectively (that is, like any MPI collective) over all processes in the graph's communicator.

Warning
It is only valid to call this overload of fillComplete if the row Map is one to one! If the row Map is NOT one to one, you must call the above three-argument version of fillComplete, and supply one-to-one domain and range Maps. If you have Maps that are not one to one, and you do not know how to make a Map that covers the same global indices but is one to one, then you may call Tpetra::createOneToOne() (see Map's header file) to make a one-to-one version of your Map.
Parameters
params[in/out] List of parameters controlling this method's behavior. See documentation of the three-argument version of fillComplete (above) for valid parameters.

Definition at line 2881 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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. 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.

This method must be called collectively (that is, like any MPI collective) over all processes in the graph's communicator.

Warning
This method is intended for expert developer use only, and should never be called by user code.
Parameters
domainMap[in] The graph's domain Map. MUST be one to one!
rangeMap[in] The graph's range Map. MUST be one to one! May be, but need not be, the same as the domain Map.
importer[in] Import from the graph's domain Map to its column Map. If no Import is necessary (i.e., if the domain and column Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
exporter[in] Export from the graph's row Map to its range Map. If no Export is necessary (i.e., if the row and range Maps are the same, in the sense of Tpetra::Map::isSameAs), then this may be Teuchos::null.
params[in/out] List of parameters controlling this method's behavior.

Definition at line 3113 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Teuchos::Comm< int > > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getComm ( ) const
overridevirtual

Returns the communicator.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 1027 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowMap ( ) const
overridevirtual

Returns the Map that describes the row distribution in this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 816 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getColMap ( ) const
overridevirtual

Returns the Map that describes the column distribution in this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 823 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getDomainMap ( ) const
overridevirtual

Returns the Map associated with the domain of this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 830 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::map_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRangeMap ( ) const
overridevirtual

Returns the Map associated with the domain of this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 837 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::import_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getImporter ( ) const
overridevirtual

Returns the importer associated with this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 844 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const typename CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::export_type > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getExporter ( ) const
overridevirtual

Returns the exporter associated with this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 851 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumRows ( ) const
overridevirtual

Returns the number of global rows in the graph.

Undefined if isFillActive().

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 776 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumCols ( ) const
overridevirtual

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 783 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalNumRows ( ) const
overridevirtual

Returns the number of graph rows owned on the calling node.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 795 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalNumCols ( ) const
overridevirtual

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 802 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
GlobalOrdinal Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getIndexBase ( ) const
overridevirtual

Returns the index base for global indices for this graph.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 1034 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalNumEntries ( ) const
overridevirtual

Returns the global number of entries in the graph.

Undefined if isFillActive ().

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 877 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalNumEntries ( ) const
overridevirtual

The local number of entries in the graph.

"Local" means "local to the calling (MPI) process."

Warning
If the graph is not fill complete, this may launch a thread-parallel computational kernel. This is because we do not store the number of entries as a separate integer field, since doing so and keeping it updated would hinder thread-parallel insertion of new entries. See #1357.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 889 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInGlobalRow ( global_ordinal_type  globalRow) const
overridevirtual

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 2018 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumEntriesInLocalRow ( local_ordinal_type  localRow) const
overridevirtual

Get the number of entries in the given row (local index).

Returns
The number of entries in the given row, specified by local index, on the calling MPI process. If the specified local row index is invalid on the calling process, return Teuchos::OrdinalTraits<size_t>::invalid().

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 2030 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalAllocationSize ( ) const

The local number of indices allocated for the graph, over all rows on the calling (MPI) process.

"Local" means "local to the calling (MPI) process."

Warning
If the graph is not fill complete, this may require computation. This is because we do not store the allocation count as a separate integer field, since doing so and keeping it updated would hinder thread-parallel insertion of new entries.

This is the allocation available to the user. Actual allocation may be larger, for example, after calling fillComplete(). Thus, this does not necessarily reflect the graph's memory consumption.

Returns
If indicesAreAllocated() is true, the allocation size. Otherwise, Tpetra::Details::OrdinalTraits<size_t>::invalid().

Definition at line 989 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumAllocatedEntriesInGlobalRow ( global_ordinal_type  globalRow) const

Current number of allocated entries in the given row on the calling (MPI) process, using a global row index.

Returns
If the given row index is in the row Map on the calling process, then return this process' allocation size for that row. Otherwise, return Tpetra::Details::OrdinalTraits<size_t>::invalid().

Definition at line 2042 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getNumAllocatedEntriesInLocalRow ( local_ordinal_type  localRow) const

Current number of allocated entries in the given row on the calling (MPI) process, using a local row index.

Returns
If the given row index is in the row Map on the calling process, then return this process' allocation size for that row. Otherwise, return Tpetra::Details::OrdinalTraits<size_t>::invalid().

Definition at line 2054 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalMaxNumRowEntries ( ) const
overridevirtual

Maximum number of entries in any row of the graph, over all processes in the graph's communicator.

Precondition
! isFillActive()
Note
This is the same as the result of a global maximum of getLocalMaxNumRowEntries() over all processes. That may not necessarily mean what you think it does if some rows of the matrix are owned by multiple processes. In particular, some processes might only own some of the entries in a particular row. This method only counts the number of entries in each row that a process owns, not the total number of entries in the row over all processes.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 946 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalMaxNumRowEntries ( ) const
overridevirtual

Maximum number of entries in any row of the graph, on this process.

Precondition
! isFillActive()

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 958 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::hasColMap ( ) const
overridevirtual

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

  • use local indices for storing entries in each row, and
  • compute an Import from the domain Map to the column Map.

The latter is mainly useful for a graph associated with a CrsMatrix.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 857 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isLocallyIndexed ( ) const
overridevirtual

Whether the graph's column indices are stored as local indices.

Weird quirk inherited from Epetra: ! isLocallyIndexed() && ! isGloballyIndexed() means that there are no graph entries on the calling process. Please don't rely on this behavior, but note that it's possible.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 976 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isGloballyIndexed ( ) const
overridevirtual

Whether the graph's column indices are stored as global indices.

Weird quirk inherited from Epetra: ! isLocallyIndexed() && ! isGloballyIndexed() means that there are no graph entries on the calling process. Please don't rely on this behavior, but note that it's possible.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 982 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isFillComplete ( ) const
overridevirtual

Whether fillComplete() has been called and the graph is in compute mode.

Implements Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 964 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isFillActive ( ) const

Whether resumeFill() has been called and the graph is in edit mode.

Definition at line 970 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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 1046 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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 863 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy ( global_ordinal_type  gblRow,
nonconst_global_inds_host_view_type &  gblColInds,
size_t &  numColInds 
) const
override

Get a copy of the given row, using global indices.

Parameters
gblRow[in] Global index of the row.
gblColInds[out] On output: Global column indices.
numColInds[out] Number of indices returned.

Definition at line 2135 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( local_ordinal_type  gblRow,
nonconst_local_inds_host_view_type &  gblColInds,
size_t &  numColInds 
) const
override

Get a copy of the given row, using local indices.

Parameters
lclRow[in] Local index of the row.
lclColInds[out] On output: Local column indices.
numColInds[out] Number of indices returned.
Precondition
hasColMap()

Definition at line 2093 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( const global_ordinal_type  gblRow,
global_inds_host_view_type &  gblColInds 
) const
override

Get a const view of the given global row's global column indices.

Parameters
gblRow[in] Global index of the row.
gblColInds[out] Global column indices in the row. If the given row is not a valid row index on the calling process, then the result has no entries (its size is zero).
Precondition
! isLocallyIndexed()
Postcondition
gblColInds.size() == getNumEntriesInGlobalRow(gblRow)

Definition at line 2204 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::supportsRowViews ( ) const
overridevirtual

Whether this class implements getLocalRowView() and getGlobalRowView() (it does).

Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 6258 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal, class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( const LocalOrdinal  lclRow,
local_inds_host_view_type &  lclColInds 
) const
override

Get a const view of the given local row's local column indices.

Parameters
lclRow[in] Local index of the row.
lclColInds[out] Local column indices in the row. If the given row is not a valid row index on the calling process, then the result has no entries (its size is zero).
Precondition
! isGloballyIndexed()
Postcondition
lclColInds.size() == getNumEntriesInLocalRow(lclRow)

Definition at line 2169 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
std::string Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::description ( ) const
overridevirtual

Return a one-line human-readable description of this object.

Reimplemented from Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 4314 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
overridevirtual

Print this object to the given output stream with the given verbosity level.

Reimplemented from Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 4333 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::checkSizes ( const SrcDistObject source)
overridevirtual

Compare the source and target (this) objects for compatibility.

Returns
True if they are compatible, else false.

Implements Tpetra::DistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 4452 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM 
)
overridevirtual

DistObject copyAndPermute has multiple overloads – use copyAndPermutes for anything we don't override

Definition at line 4461 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::pack ( const Teuchos::ArrayView< const local_ordinal_type > &  exportLIDs,
Teuchos::Array< global_ordinal_type > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets 
) const
overridevirtual

< DistObject overloads packAndPrepare. Explicitly use DistObject's packAndPrepare for anything we don't override

Reimplemented from Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 5233 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  importLIDs,
Kokkos::DualView< packet_type *, buffer_device_type imports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
const size_t  constantNumPackets,
const CombineMode  combineMode 
)
overridevirtual

< DistObject has overloaded unpackAndCombine, use the DistObject's implementation for anything we don't override.

Definition at line 5660 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagOffsets ( const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &  offsets) const

Get offsets of the diagonal entries in the graph.

Warning
This method is only for expert users.
We make no promises about backwards compatibility for this method. It may disappear or change at any time.
This method must be called collectively. We reserve the right to do extra checking in a debug build that will require collectives.

This method helps users optimize Tpetra::CrsMatrix::getLocalDiagCopy and Tpetra::Experimental::BlockCrsMatrix::getLocalDiagCopy, for several calls when the graph's structure does not change. The method fills an array of offsets of the local diagonal entries in the matrix. getLocalDiagCopy uses the offsets to extract the diagonal entries directly, without needing to search for them using Map lookups and search in each row of the graph.

The output array's contents are not defined in any other context other than for use in getLocalDiagCopy. For example, you should not rely on offsets(i) being the index of the diagonal entry in the views returned by Tpetra::CrsMatrix::getLocalRowView. This may be the case, but it need not be. (For example, we may choose to optimize the lookups down to the optimized storage level, in which case the offsets will be computed with respect to the underlying storage format, rather than with respect to the views.)

Changes to the graph's structure, or calling fillComplete on the graph (if its structure is not already fixed), may make the output array's contents invalid. "Invalid" means that you must call this method again to recompute the offsets.

Precondition
The graph must have a column Map.
All diagonal entries of the graph must be populated on this process. Results are undefined otherwise.
offsets.extent(0) >= this->getLocalNumRows()
Parameters
offsets[out] Output array of offsets. This method does NOT allocate the array; the caller must allocate. Must have getLocalNumRows() entries on the calling process. (This may be different on different processes.)

Definition at line 5910 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalOffRankOffsets ( offset_device_view_type &  offsets) const

Get offsets of the off-rank entries in the graph.

Definition at line 6096 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalDiagOffsets ( Teuchos::ArrayRCP< size_t > &  offsets) const

Backwards compatibility overload of the above method.

This method takes a Teuchos::ArrayRCP instead of a Kokkos::View. It also reallocates the output array if it is not long enough.

Parameters
offsets[out] Output array of offsets. This method reallocates the array if it is not long enough. This is why the method takes the array by reference.

Definition at line 6223 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setAllIndices ( const typename local_graph_device_type::row_map_type &  rowPointers,
const typename local_graph_device_type::entries_type::non_const_type &  columnIndices 
)

Set the graph's data directly, using 1-D storage.

Precondition
columnIndices are sorted within rows
hasColMap() == true
rowPointers.size() == getLocalNumRows()+1
No insert routines have been called.
Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 2496 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::setAllIndices ( const Teuchos::ArrayRCP< size_t > &  rowPointers,
const Teuchos::ArrayRCP< local_ordinal_type > &  columnIndices 
)

Set the graph's data directly, using 1-D storage.

Precondition
columnIndices are sorted within rows
hasColMap() == true
rowPointers.size() == getLocalNumRows()+1
No insert routines have been called.
Warning
This method is intended for expert developer use only, and should never be called by user code.

Definition at line 2602 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::row_ptrs_host_view_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowPtrsHost ( ) const

Get a host view of the packed row offsets.

Definition at line 2066 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::row_ptrs_device_view_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowPtrsDevice ( ) const

Get a device view of the packed row offsets.

Definition at line 2073 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_host_view_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalIndicesHost ( ) const

Get a host view of the packed column indicies.

Definition at line 2080 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_device_view_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalIndicesDevice ( ) const

Get a device view of the packed column indicies.

Definition at line 2087 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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.

Parameters
newColMap[in] New column Map. Must be nonnull. Within Tpetra, there are no particular restrictions on the column map. However, if this graph will be used in Xpetra, Ifpack2, or MueLu, the column map's list of global indices must follow "Aztec ordering": locally owned GIDs (same order as in domain map), followed by remote GIDs (in order of owning proc, and sorted within each proc).

It is strongly recommended to use Tpetra::Details::makeColMap() to create the column map. makeColMap() follows Aztec ordering by default.

Definition at line 3499 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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.

Precondition
On every calling process, every index owned by the current column Map must also be owned by the new column Map.
If the new Import object is provided, the new Import object's source Map must be the same as the current domain Map, and the new Import's target Map must be the same as the new column Map.
Parameters
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 3516 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::replaceDomainMap ( const Teuchos::RCP< const map_type > &  newDomainMap)

Replace the current domain Map with the given objects.

The Graph's Import object will be recomputed if needed.

Precondition
isFillComplete() == true
isFillActive() == false

Definition at line 3746 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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.

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

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::replaceRangeMap ( const Teuchos::RCP< const map_type > &  newRangeMap)

Replace the current Range Map with the given objects.

The Graph's Export object will be recomputed if needed.

Precondition
isFillComplete() == true
isFillActive() == false

Definition at line 3804 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::replaceRangeMapAndExporter ( const Teuchos::RCP< const map_type > &  newRangeMap,
const Teuchos::RCP< const export_type > &  newExporter 
)

Replace the current Range Map and Export with the given parameters.

Warning
This method is ONLY for use by experts.
We make NO promises of backwards compatibility. This method may change or disappear at any time.
Precondition
isFillComplete() == true
isFillActive() == false
Either the given Export object is null, or the target Map of the given Export is the same as this graph's Range Map.
Either the given Export object is null, or the source Map of the given Export is the same as this graph's Row Map.

Definition at line 3823 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::removeEmptyProcessesInPlace ( const Teuchos::RCP< const map_type > &  newMap)
overridevirtual

Remove processes owning zero rows from the Maps and their communicator.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in Tpetra_DistObject_decl.hpp.
We make NO promises of backwards compatibility. This method may change or disappear at any time.
Parameters
newMap[in] This must be the result of calling the removeEmptyProcesses() method on the row Map. If it is not, this method's behavior is undefined. This pointer will be null on excluded processes.

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 5815 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node>
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::importAndFillComplete ( Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &  destGraph,
const import_type importer,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
) const

Import from this to the given destination graph, and make the result fill complete.

If destGraph.is_null(), this creates a new graph as the destination. (This is why destGraph is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the graph has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsGraph, is preferred for user applications.

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

Definition at line 6866 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node>
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::importAndFillComplete ( Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &  destGraph,
const import_type rowImporter,
const import_type domainImporter,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
) const

Import from this to the given destination graph, and make the result fill complete.

If destGraph.is_null(), this creates a new graph as the destination. (This is why destGraph is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the graph has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsGraph, is preferred for user applications.

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

Definition at line 6876 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node>
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::exportAndFillComplete ( Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &  destGraph,
const export_type exporter,
const Teuchos::RCP< const map_type > &  domainMap = Teuchos::null,
const Teuchos::RCP< const map_type > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
) const

Export from this to the given destination graph, and make the result fill complete.

If destGraph.is_null(), this creates a new graph as the destination. (This is why destGraph is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the graph has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsGraph, is preferred for user applications.

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

Definition at line 6887 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node>
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::exportAndFillComplete ( Teuchos::RCP< CrsGraph< local_ordinal_type, global_ordinal_type, Node >> &  destGraph,
const export_type rowExporter,
const export_type domainExporter,
const Teuchos::RCP< const map_type > &  domainMap,
const Teuchos::RCP< const map_type > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
) const

Export from this to the given destination graph, and make the result fill complete.

If destGraph.is_null(), this creates a new graph as the destination. (This is why destGraph is passed in by nonconst reference to RCP.) Otherwise it checks for "pristine" status and throws if that is not the case. "Pristine" means that the graph has no entries and is not fill complete.

Use of the "non-member constructor" version of this method, exportAndFillCompleteCrsGraph, is preferred for user applications.

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

Definition at line 6897 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::makeColMap ( Teuchos::Array< int > &  remotePIDs)
protected

Make and set the graph's column Map.

This method makes the column Map, even if the graph already has one. It is the caller's responsibility not to call this method unnecessarily.

Parameters
remotePIDs[out] The process ranks corresponding to the column Map's "remote" (not on the calling process in the domain Map) indices.

Definition at line 4106 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
std::pair< size_t, std::string > Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::makeIndicesLocal ( const bool  verbose = false)
protected

Convert column indices from global to local.

Precondition
The graph has a column Map.
Postcondition
The graph is locally indexed.
Parameters
verbose[in] Whether to print verbose debugging output. This exists because CrsMatrix may want to control output independently of the CrsGraph that it owns.
Returns
Error code and error string. See below.

First return value is the number of column indices on this process, counting duplicates, that could not be converted to local indices, because they were not in the column Map on the calling process. If some error occurred before conversion happened, then this is Tpetra::Details::OrdinalTraits<size_t>::invalid().

Second return value is a human-readable error string. If the first return value is zero, then the string may be empty.

Definition at line 3950 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::makeImportExport ( Teuchos::Array< int > &  remotePIDs,
const bool  useRemotePIDs 
)
protected

Make the Import and Export objects, if needed.

Parameters
remotePIDs[in/out] On input: the output of makeColMap(). May be modified on output.
useRemotePIDs[in] Whether to use remotePIDs. Use it if we called makeColMap with this as the output argument, else don't use it.

Definition at line 4253 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertIndices ( RowInfo rowInfo,
const SLocalGlobalViews &  newInds,
const ELocalGlobal  lg,
const ELocalGlobal  I 
)
protected

Insert indices into the given row.

Precondition
! (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.

This method does no allocation; it just inserts the indices.

Parameters
rowInfo[in/out] On input: 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. On output: numEntries field is updated.
newInds[in] View of the column indices to insert. If lg == GlobalIndices, then newInds.ginds, a Teuchos::ArrayView<const global_ordinal_type>, contains the (global) column indices to insert. Otherwise, if lg == LocalIndices, then newInds.linds, a Teuchos::ArrayView<const local_ordinal_type>, contains the (local) column indices to insert.
lgIf lg == GlobalIndices, then the input indices (in newInds) are global indices. Otherwise, if lg == LocalIndices, the input indices are local indices.
IIf 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.
Returns
The number of indices inserted.

Definition at line 1448 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndicesImpl ( const local_ordinal_type  lclRow,
const global_ordinal_type  inputGblColInds[],
const size_t  numInputInds 
)
protected

Insert global indices, using an input local row index.

Parameters
rowInfo[in] Result of getRowInfo() on the row in which to insert.
inputGblColInds[in] Input global column indices.
numInputInds[in] The number of input global column indices.
Returns
The number of indices inserted.

Definition at line 1545 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndicesImpl ( const RowInfo rowInfo,
const global_ordinal_type  inputGblColInds[],
const size_t  numInputInds,
std::function< void(const size_t, const size_t, const size_t)>  fun = std::function<void(const size_t, const size_t, const size_t)>() 
)
protected

Insert global indices, using an input RowInfo.

Parameters
rowInfo[in] Result of getRowInfo() on the row in which to insert.
inputGblColInds[in] Input global column indices.
numInputInds[in] The number of input global column indices.
Returns
The number of indices inserted.

Definition at line 1555 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::findGlobalIndices ( const RowInfo rowInfo,
const Teuchos::ArrayView< const global_ordinal_type > &  indices,
std::function< void(const size_t, const size_t, const size_t)>  fun 
) const
protected

Finds indices in the given row.

This method does no insertion; it just finds indices and calls a callback for each found index

Parameters
row[in] Row of interest
indices[in] Column indices to find in row
funCall back function called at each found index. Called as fun(k, start, offset); where k is the index in to indices, start is offset to the start of the row, and offset is the relative offset of indices[k] in the graphs indices.
Returns
The number of indices found.

Definition at line 1682 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndicesFiltered ( const local_ordinal_type  lclRow,
const global_ordinal_type  gblColInds[],
const local_ordinal_type  numGblColInds 
)
protected

Like insertGlobalIndices(), but with column Map filtering.

"Column Map filtering" means that any column indices not in the column Map on the calling process, get silently dropped.

Parameters
lclRow[in] Local index of the row in which to insert. This row MUST be in the row Map on the calling process.
gblColInds[in] The global column indices to insert into that row.
numGblColInds[in] The number of global column indices to insert into that row.

Definition at line 2389 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::insertGlobalIndicesIntoNonownedRows ( const global_ordinal_type  gblRow,
const global_ordinal_type  gblColInds[],
const local_ordinal_type  numGblColInds 
)
protected

Implementation of insertGlobalIndices for nonowned rows.

A global row index is nonowned when it is not in the column Map on the calling process.

Parameters
gblRow[in] Global index of the row in which to insert. This row must NOT be in the row Map on the calling process.
gblColInds[in] The global column indices to insert into that row.
numGblColInds[in] The number of global column indices to insert into that row.

Definition at line 2448 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::isMerged ( ) const
protected

Whether duplicate column indices in each row have been merged.

Definition at line 1052 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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 1058 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortAndMergeRowIndices ( const RowInfo rowInfo,
const bool  sorted,
const bool  merged 
)
protected

Sort and merge duplicate column indices in the given row.

Precondition
The graph is locally indexed: isGloballyIndexed() == false.
The graph is not already storage optimized: isStorageOptimized() == false
Returns
The number of duplicate column indices eliminated from the row.

Definition at line 1714 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
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).

Parameters
domainMap[in] The new domain Map
rangeMap[in] The new range Map

Definition at line 1746 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants ( ) const
inline

Returns true if globalConstants have been computed; false otherwise.

Definition at line 2008 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::computeGlobalConstants ( )

Compute global constants, if they have not yet been computed.

Warning
This is an implementation detail of Tpetra. It may change or disappear at any time. It is public only because MueLu setup needs it to be public.

Global constants include:

  • globalNumEntries_
  • globalMaxNumRowEntries_

Always compute the following:

  • globalNumEntries_
  • globalMaxNumRowEntries_

Definition at line 3880 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::computeLocalConstants ( )
protected

Compute local constants, if they have not yet been computed.

Warning
You MUST call fillLocalGraph (or CrsMatrix::fillLocalGraphAndMatrix) before calling this method! This method depends on the Kokkos::StaticCrsGraph (local_graph_device_type) object being ready.

Local constants include:

  • nodeMaxNumRowEntries_

Always compute the following:

  • nodeMaxNumRowEntries_

computeGlobalConstants calls this method, if global constants have not yet been computed.

Definition at line 3923 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
RowInfo Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowInfo ( const local_ordinal_type  myRow) const
protected

Get information about the locally owned row with local index myRow.

Definition at line 1315 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
RowInfo Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowInfoFromGlobalRowIndex ( const global_ordinal_type  gblRow) const
protected

Get information about the locally owned row with global index gblRow.

If gblRow is not locally owned, the localRow field of the returned struct is Teuchos::OrdinalTraits<size_t>::invalid().

The point of this method is to fuse the global-to-local row index lookup for checking whether gblRow is locally owned, with other error checking that getRowInfo() does. This avoids an extra global-to-local index lookup in methods like CrsMatrix::replaceGlobalValues().

Definition at line 1355 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_graph_device_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalGraphDevice ( ) const

Get the local graph.

Warning
THIS IS AN EXPERT MODE FUNCTION. THIS IS AN IMPLEMENTATION DETAIL. DO NOT CALL THIS FUNCTION!!!

This is only a valid representation of the local graph if the (global) graph is fill complete.

Definition at line 3863 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::checkInternalState ( ) const
protected

Throw an exception if the internal state is not consistent.

Definition at line 1771 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node>
void Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::swap ( CrsGraph< local_ordinal_type, global_ordinal_type, Node > &  graph)
protected

Swaps the data from *this with the data and maps from graph.

Parameters
graph[in/out] a crsGraph

Definition at line 6908 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
const row_ptrs_device_view_type& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowPtrsUnpackedDevice ( ) const
inlineprotected

Get the unpacked row pointers on device.

Definition at line 2173 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
const row_ptrs_host_view_type& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowPtrsUnpackedHost ( ) const
inlineprotected

Get the unpacked row pointers on host. Lazily make a copy from device.

Definition at line 2178 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
const row_ptrs_device_view_type& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowPtrsPackedDevice ( ) const
inlineprotected

Get the packed row pointers on device.

Definition at line 2209 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
const row_ptrs_host_view_type& Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getRowPtrsPackedHost ( ) const
inlineprotected

Get the packed row pointers on host. Lazily make a copy from device.

Definition at line 2214 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_dualv_type::t_host::const_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalIndsViewHost ( const RowInfo rowinfo) const
protected

Get a const, locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 1251 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_dualv_type::t_dev::const_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalIndsViewDevice ( const RowInfo rowinfo) const
protected

Get a const, locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 1290 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::global_inds_dualv_type::t_host::const_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalIndsViewHost ( const RowInfo rowinfo) const
protected

Get a const, globally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 1277 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::global_inds_dualv_type::t_dev::const_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalIndsViewDevice ( const RowInfo rowinfo) const
protected

Get a const, globally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 1303 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::local_inds_dualv_type::t_host Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalIndsViewHostNonConst ( const RowInfo rowinfo)
protected

Get a ReadWrite locally indexed view of the locally owned row myRow, such that rowinfo = getRowInfo(myRow).

Definition at line 1264 of file Tpetra_CrsGraph_def.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowCopy ( const GlobalOrdinal  gblRow,
nonconst_global_inds_host_view_type &  gblColInds,
size_t &  numColInds 
) const
pure virtualinherited

Get a copy of the global column indices in a given row of the graph.

Given the global index of a row of the graph, get a copy of all the global column indices in that row that the calling process stores.

Parameters
gblRow[in] Global index of the row.
gblColInds[in/out] On output: All the global column indices in that row on the calling process.
numColInds[out] Number of indices in the row on the calling process.
Precondition
getRowMap()->isNodeGlobalElement(gblRow) is true.
gblColInds.size() >= getNumEntriesInGlobalRow(gblRow) is true.
template<class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowCopy ( const LocalOrdinal  lclRow,
nonconst_local_inds_host_view_type &  lclColInds,
size_t &  numColInds 
) const
pure virtualinherited

Get a copy of the local column indices in a given row of the graph.

Given the local index of a row of the graph, get a copy of all the local column indices in that row that the calling process stores.

Parameters
lclRow[in] Local index of the row.
lclColInds[in/out] On output: All the local column indices in that row on the calling process.
numColInds[out] Number of indices in the row on the calling process.
Precondition
hasColMap() is true.
getRowMap()->isNodeLocalElement(lclRow) is true.
lclColInds.size() >= getNumEntriesInLocalRow(lclRow) is true.
template<class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >::getLocalRowView ( const LocalOrdinal  lclRow,
local_inds_host_view_type &  lclColInds 
) const
pure virtualinherited

Get a constant, nonpersisting, locally indexed view of the given row of the graph.

The returned views of the column indices are not guaranteed to persist beyond the lifetime of this. Furthermore, some RowGraph implementations allow changing the values, or the indices and values. Any such changes invalidate the returned views.

This method only gets the entries in the given row that are stored on the calling process. Note that if the graph has an overlapping row Map, it is possible that the calling process does not store all the entries in that row.

Precondition
isLocallyIndexed () && supportsRowViews ()
Postcondition
indices.size () == getNumEntriesInGlobalRow (LocalRow)
Parameters
lclRow[in] Local index of the row.
lclColInds[out] Local indices of the columns in the row. If the given row is not a valid row index on the calling process, then the result has no entries (its size is zero).

Subclasses are expected to implement this method. We would have made this method pure virtual, but that would have broken backwards compatibility, since we added the method at least one major release after introducing this class.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
virtual void Tpetra::RowGraph< LocalOrdinal, GlobalOrdinal, Node >::getGlobalRowView ( const GlobalOrdinal  gblRow,
global_inds_host_view_type &  gblColInds 
) const
pure virtualinherited

Get a const, non-persisting view of the given global row's global column indices, as a Teuchos::ArrayView.

Parameters
gblRow[in] Global index of the row.
gblColInds[out] Global column indices in the row. If the given row is not a valid row index on the calling process, then the result has no entries (its size is zero).
Precondition
! isLocallyIndexed()
Postcondition
gblColInds.size() == getNumEntriesInGlobalRow(gblRow)

Subclasses are expected to implement this method. We would have made this method pure virtual, but that would have broken backwards compatibility, since we added the method at least one major release after introducing this class.

void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::doImport ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
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,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
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,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
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,
const CombineMode  CM,
const bool  restrictedMode = false 
)
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.

"Restricted Mode" does two things:

  1. Skips copyAndPermute
  2. Allows the "target" Map of the transfer to be a subset of the Map of *this, in a "locallyFitted" sense.

This cannot be used if (2) is not true, OR there are permutes. The "source" maps still need to match.

Parameters
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 >::transferArrived ( ) const
inherited

Whether the data from an import/export operation has arrived, and is ready for the unpack and combine step.

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_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::getMap ( ) const
inlinevirtualinherited

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 554 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_type > &  newMap)
virtualinherited

Remove processes which contain no entries in this object's Map.

Warning
This method is ONLY for use by experts. We highly recommend using the nonmember function of the same name defined in this file.
We make NO promises of backwards compatibility. This method may change or disappear at any time.

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 entries ("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.

Returns
The object's new Map. Its communicator is a new communicator, distinct from the old Map's communicator, which contains a subset of the processes in the old communicator.
Note
The name differs from Map's method removeEmptyProcesses(), in order to emphasize that the operation on DistObject happens in place, modifying the input, whereas the operation removeEmptyProcess() on Map does not modify the input.
To implementers of DistObject subclasses: The default implementation of this class throws std::logic_error.
const Details::DistributorActor& Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::getActor ( )
inlineinherited

Getter for the DistributorActor.

Definition at line 639 of file Tpetra_DistObject_decl.hpp.

virtual size_t Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::constantNumberOfPackets ( ) const
protectedvirtualinherited

Whether the implementation's instance promises always to have a constant number of packets per LID (local index), 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 (local index) 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,
const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &  transfer,
const char  modeString[],
const ReverseOption  revOp,
const CombineMode  CM,
const bool  restrictedMode 
)
protectedvirtualinherited

Redistribute data across (MPI) processes.

Parameters
src[in] The source object, to redistribute into the target object, which is *this object.
transfer[in] The Export or Import object representing the communication pattern. (Details::Transfer is the common base class of these two objects.)
modeString[in] Human-readable string, for verbose debugging output and error output, explaining what function called this method. Example: "doImport (forward)", "doExport (reverse)".
revOp[in] Whether to do a forward or reverse mode redistribution.
CM[in] The combine mode that describes how to combine values that map to the same global ID on the same process.
virtual bool Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::reallocArraysForNumPacketsPerLid ( const size_t  numExportLIDs,
const size_t  numImportLIDs 
)
protectedvirtualinherited

Reallocate numExportPacketsPerLID_ and/or numImportPacketsPerLID_, if necessary.

Parameters
numExportLIDs[in] Number of entries in the exportLIDs input array argument of doTransfer().
numImportLIDs[in] Number of entries in the remoteLIDs input array argument of doTransfer().
Returns
Whether we actually reallocated either of the arrays.
Warning
This is an implementation detail of doTransferNew(). This needs to be protected, but that doesn't mean users should call this method.
void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::beginTransfer ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  src,
const ::Tpetra::Details::Transfer< local_ordinal_type, global_ordinal_type, node_type > &  transfer,
const char  modeString[],
const ReverseOption  revOp,
const CombineMode  CM,
const bool  restrictedMode 
)
protectedinherited

Implementation detail of doTransfer.

LID DualViews come from the Transfer object given to doTransfer. They are always sync'd on both host and device. Users must never attempt to modify or sync them.

virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM 
)
protectedvirtualinherited

Perform copies and permutations that are local to the calling (MPI) process.

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, packs the source object's data.

Precondition
permuteToLIDs and permuteFromLIDs are sync'd to both host and device. That is, permuteToLIDs.need_sync_host(), permuteToLIDs.need_sync_device(), permuteFromLIDs.need_sync_host(), and permuteFromLIDs.need_sync_device() are all false.
Parameters
source[in] On entry, the source object of the Export or Import operation.
numSameIDs[in] The number of elements that are the same on the source and target objects. These elements live on the same process in both the source and target objects.
permuteToLIDs[in] List of the elements that are permuted. They are listed by their local index (LID) in the destination object.
permuteFromLIDs[in] List of the elements that are permuted. They are listed by their local index (LID) in the source object.
CM[in] CombineMode to be used during copyAndPermute; may or may not be used by the particular object being called; behavior with respect to CombineMode may differ by object.
virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::copyAndPermute ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs,
const CombineMode  CM,
const execution_space space 
)
protectedvirtualinherited

Same as copyAndPermute, but do operations in space.

virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< packet_type *, buffer_device_type > &  exports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
size_t &  constantNumPackets 
)
protectedvirtualinherited

Pack data and metadata for communication (sends).

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, packs the source object's data.

Precondition
exportLIDs is sync'd to both host and device. That is, exportLIDs.need_sync_host () and exportLIDs.need_sync_device() are both false.
Parameters
source[in] Source object for the redistribution.
exportLIDs[in] List of the entries (as local IDs in the source object) that Tpetra will send to other processes.
exports[out] On exit, the packed data to send. Implementations must reallocate this as needed (prefer reusing the existing allocation if possible), and may modify and/or sync this wherever they like.
numPacketsPerLID[out] On exit, the implementation of this method must do one of two things: either set numPacketsPerLID[i] to the number of packets to be packed for exportLIDs[i] and set constantNumPackets to zero, or set constantNumPackets to a nonzero value. If the latter, the implementation must not modify the entries of numPacketsPerLID. If the former, the implementation may sync numPacketsPerLID this wherever it likes, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
constantNumPackets[out] On exit, 0 if the number of packets per LID could differ, else (if nonzero) the number of packets per LID (which must be constant).
virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::packAndPrepare ( const SrcDistObject< GlobalOrdinal, LocalOrdinal, GlobalOrdinal, Node > &  source,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< packet_type *, buffer_device_type > &  exports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
size_t &  constantNumPackets,
const execution_space space 
)
protectedvirtualinherited

Same as packAndPrepare, but in an execution space instance.

virtual void Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::unpackAndCombine ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  importLIDs,
Kokkos::DualView< packet_type *, buffer_device_type imports,
Kokkos::DualView< size_t *, buffer_device_type numPacketsPerLID,
const size_t  constantNumPackets,
const CombineMode  combineMode 
)
protectedvirtualinherited

Perform any unpacking and combining after communication.

Subclasses must reimplement this function. Its default implementation does nothing. Note that the <t>target object of the Export or Import, namely *this, unpacks the received data into itself, possibly modifying its entries.

Precondition
importLIDs is sync'd to both host and device. That is, importLIDs.need_sync_host () and importLIDs.need_sync_device() are both false.
Parameters
importLIDs[in] List of the entries (as LIDs in the destination object) we received from other processes.
imports[in/out] On input: Buffer of received data to unpack. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
numPacketsPerLID[in/out] On input: If constantNumPackets is zero, then numPacketsPerLID[i] contains the number of packets imported for importLIDs[i]. DistObject promises nothing about where this is sync'd. Implementations may sync this wherever they like, either to host or to device. The allocation belongs to DistObject, not to subclasses; don't be tempted to change this to pass by reference.
constantNumPackets[in] If nonzero, then the number of packets per LID is the same for all entries ("constant") and constantNumPackets is that number. If zero, then numPacketsPerLID[i] is the number of packets to unpack for LID importLIDs[i].
combineMode[in] The CombineMode to use when combining the imported entries with existing entries.
virtual bool Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::reallocImportsIfNeeded ( const size_t  newSize,
const bool  verbose,
const std::string *  prefix,
const bool  remoteLIDsContiguous = false,
const CombineMode  CM = INSERT 
)
protectedvirtualinherited

Reallocate imports_ if needed.

This unfortunately must be declared protected, for the same reason that imports_ is declared protected.

Parameters
newSize[in] New size of imports_.
verbose[in] Whether to print verbose debugging output to stderr on every (MPI) process in the communicator.
prefix[in] If verbose is true, then this is a nonnull prefix to print at the beginning of each line of verbose debugging output. Otherwise, not used.
Returns
Whether we actually reallocated.

We don't need a "reallocExportsIfNeeded" method, because exports_ always gets passed into packAndPrepare() by nonconst reference. Thus, that method can resize the DualView without needing to call other DistObject methods.

Friends And Related Function Documentation

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 > &  params = Teuchos::null 
)
related

Nonmember function to create an empty CrsGraph given a row Map and the max number of entries allowed locally per row.

Returns
A graph with at most the specified number of nonzeros per row (defaults to zero).

Definition at line 2524 of file Tpetra_CrsGraph_decl.hpp.

template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph ( const Teuchos::RCP< const CrsGraphType > &  sourceGraph,
const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  importer,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
related

Nonmember CrsGraph constructor that fuses Import and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

A common use case is to create an empty destination CrsGraph, redistribute from a source CrsGraph (by an Import or Export operation), then call fillComplete() on the destination CrsGraph. This constructor fuses these three cases, for an Import redistribution.

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsGraph.

The resulting graph is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source graph, and its range Map is the range Map of the source graph.

Warning
If the target Map of the Import is a subset of the source Map of the Import, then you cannot use the default range Map. You should instead construct a nonoverlapping version of the target Map and supply that as the nondefault value of the range Map.
Parameters
sourceGraph[in] The source graph from which to import. The source of an Import must have a nonoverlapping distribution.
importer[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the rowMap of sourceGraph unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the rowMap of the sourceGraph
domainMap[in] Domain Map of the returned graph. If null, we use the default, which is the domain Map of the source graph.
rangeMap[in] Range Map of the returned graph. If null, we use the default, which is the range Map of the source graph.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 2587 of file Tpetra_CrsGraph_decl.hpp.

template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph ( const Teuchos::RCP< const CrsGraphType > &  sourceGraph,
const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  rowImporter,
const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  domainImporter,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  domainMap,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)
related

Nonmember CrsGraph constructor that fuses Import and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

A common use case is to create an empty destination CrsGraph, redistribute from a source CrsGraph (by an Import or Export operation), then call fillComplete() on the destination CrsGraph. This constructor fuses these three cases, for an Import redistribution.

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsGraph.

The resulting graph is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source graph, and its range Map is the range Map of the source graph.

Warning
If the target Map of the Import is a subset of the source Map of the Import, then you cannot use the default range Map. You should instead construct a nonoverlapping version of the target Map and supply that as the nondefault value of the range Map.
Parameters
sourceGraph[in] The source graph from which to import. The source of an Import must have a nonoverlapping distribution.
rowImporter[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the rowMap of sourceGraph unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the rowMap of the sourceGraph
domainImporter[in] The Import instance containing a precomputed redistribution plan. The source Map of the Import must be the same as the domainMap of sourceGraph unless the "Reverse Mode" option on the params list, in which case the targetMap of Import must match the domainMap of the sourceGraph
domainMap[in] Domain Map of the returned graph.
rangeMap[in] Range Map of the returned graph.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 2655 of file Tpetra_CrsGraph_decl.hpp.

template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph ( const Teuchos::RCP< const CrsGraphType > &  sourceGraph,
const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  exporter,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)
related

Nonmember CrsGraph constructor that fuses Export and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

For justification, see the documentation of importAndFillCompleteCrsGraph() (which is the Import analog of this function).

The resulting graph is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source graph, and its range Map is the range Map of the source graph.

Parameters
sourceGraph[in] The source graph from which to export. Its row Map may be overlapping, since the source of an Export may be overlapping.
exporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the row Map of sourceGraph.
domainMap[in] Domain Map of the returned graph. If null, we use the default, which is the domain Map of the source graph.
rangeMap[in] Range Map of the returned graph. If null, we use the default, which is the range Map of the source graph.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 2709 of file Tpetra_CrsGraph_decl.hpp.

template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph ( const Teuchos::RCP< const CrsGraphType > &  sourceGraph,
const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  rowExporter,
const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &  domainExporter,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  domainMap,
const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type >> &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)
related

Nonmember CrsGraph constructor that fuses Export and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

For justification, see the documentation of importAndFillCompleteCrsGraph() (which is the Import analog of this function).

The resulting graph is fill complete (in the sense of isFillComplete()) and has optimized storage (in the sense of isStorageOptimized()). By default, its domain Map is the domain Map of the source graph, and its range Map is the range Map of the source graph.

Parameters
sourceGraph[in] The source graph from which to export. Its row Map may be overlapping, since the source of an Export may be overlapping.
rowExporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the row Map of sourceGraph.
domainExporter[in] The Export instance containing a precomputed redistribution plan. The source Map of the Export must be the same as the domain Map of sourceGraph.
domainMap[in] Domain Map of the returned graph.
rangeMap[in] Range Map of the returned graph.
params[in/out] Optional list of parameters. If not null, any missing parameters will be filled in with their default values.

Definition at line 2760 of file Tpetra_CrsGraph_decl.hpp.

Member Data Documentation

template<class LocalOrdinal, class GlobalOrdinal, class Node>
const bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::useAtomicUpdatesByDefault
staticprotected
Initial value:
=
true

Whether transformLocalValues should use atomic updates by default.

Warning
This is an implementation detail.

Definition at line 1945 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
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 2097 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
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 2099 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
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 2101 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
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 2103 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP<const import_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::importer_
protected

The Import from the domain Map to the column Map.

This gets constructed by fillComplete. It may be null if the domain Map and the column Map are the same, since no Import is necessary in that case for sparse matrix-vector multiply.

Definition at line 2111 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Teuchos::RCP<const export_type> Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::exporter_
protected

The Export from the row Map to the range Map.

This gets constructed by fillComplete. It may be null if the row Map and the range Map are the same, since no Export is necessary in that case for sparse matrix-vector multiply.

Definition at line 2118 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::nodeMaxNumRowEntries_
protected
Initial value:
=
Teuchos::OrdinalTraits<size_t>::invalid()

Local maximum of the number of entries in each row.

Computed in computeLocalConstants; only valid when isFillComplete() is true.

Definition at line 2124 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::globalNumEntries_
protected
Initial value:
=
Teuchos::OrdinalTraits<global_size_t>::invalid()

Global number of entries in the graph.

Only valid when isFillComplete() is true.

Definition at line 2130 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
global_size_t Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::globalMaxNumRowEntries_
protected
Initial value:
=
Teuchos::OrdinalTraits<global_size_t>::invalid()

Global maximum of the number of entries in each row.

Computed in computeGlobalConstants(); only valid when isFillComplete() is true.

Definition at line 2137 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
local_inds_wdv_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::lclIndsUnpacked_wdv
protected

Local ordinals of column indices for all rows Valid when isLocallyIndexed is true If OptimizedStorage, storage is PACKED after fillComplete If not OptimizedStorate, storage is UNPACKED after fillComplete; that is, the views have storage equal to sizes provided in CrsGraph constructor.

This is allocated only if

  • The calling process has a nonzero number of entries
  • The graph is locally indexed

UVM Removal Note: Device view takes place of k_lclInds1D_

Definition at line 2265 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
local_inds_wdv_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::lclIndsPacked_wdv
mutableprotected

Local ordinals of column indices for all rows Valid when isLocallyIndexed is true Built during fillComplete or non-fillComplete constructors Storage is PACKED after fillComplete that is, the views have storage equal to sizes provided in CrsGraph constructor.

This is allocated only if

  • The calling process has a nonzero number of entries
  • The graph is locally indexed

UVM Removal Note: Device view takes place of lclGraph_.entries

Definition at line 2280 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
global_inds_wdv_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::gblInds_wdv
protected

Global ordinals of column indices for all rows.

Global ordinals of column indices for all rows

This is allocated only if

  • The calling process has a nonzero number of entries
  • The graph is globally indexed

UVM Removal Note: Device view takes place of k_gblInds1D_

Definition at line 2295 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Kokkos::View<const size_t*, device_type>::host_mirror_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::k_numAllocPerRow_
protected

The maximum number of entries to allow in each locally owned row, per row.

This comes in as an argument to some of the graph's constructors. Either this or numAllocForAllRows_ is used, but not both. allocateIndices(), setAllIndices(), and expertStaticFillComplete() all deallocate this array once they are done with it.

This is a host View because it is only ever used on the host. It has the host_mirror_type type for backwards compatibility; this used to be a DualView with default layout, so making this a host_mirror_type ensures that we can still take it directly by assignment from the constructors that take DualView, without a deep copy.

This array only exists on a process before the graph's indices are allocated on that process. After that point, it is discarded, since the graph's allocation implicitly or explicitly represents the same information.

FIXME (mfh 07 Aug 2014) We want graph's constructors to allocate, rather than doing lazy allocation at first insert. This will make both k_numAllocPerRow_ and numAllocForAllRows_ obsolete.

Definition at line 2357 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
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 k_numAllocPerRow_ is used, but not both.

FIXME (mfh 07 Aug 2014) We want graph's constructors to allocate, rather than doing lazy allocation at first insert. This will make both k_numAllocPerRow_ and numAllocForAllRows_ obsolete.

Definition at line 2368 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
num_row_entries_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::k_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.

This may also exist with 1-D storage, if storage is unpacked.

Definition at line 2420 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
offset_device_view_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::k_offRankOffsets_
mutableprotected

The offsets for off-rank entries.

When off-rank entries are sorted last, this rowPtr-lile view contains the offsets. It is compute on the first call to getLocalOffRankOffsets().

Definition at line 2427 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Details::EStorageStatus Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::storageStatus_
protected
Initial value:
=
Details::STORAGE_1D_UNPACKED

Status of the graph's storage, when not in a fill-complete state.

The phrase "When not in a fill-complete state" is important. When the graph is fill complete, it always uses 1-D "packed" storage. However, if the "Optimize Storage" parameter to fillComplete was false, the graph may keep unpacked 1-D storage around and resume it on the next resumeFill call.

Definition at line 2440 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::indicesAreSorted_ = true
protected

Whether the graph's indices are sorted in each row, on this process.

Definition at line 2449 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::noRedundancies_ = true
protected

Whether the graph's indices are non-redundant (merged) in each row, on this process.

Definition at line 2452 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::haveLocalConstants_ = false
protected

Whether this process has computed local constants.

Definition at line 2454 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::haveGlobalConstants_ = false
protected

Whether all processes have computed global constants.

Definition at line 2456 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
nonlocals_type Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::nonlocals_
protected

Nonlocal data given to insertGlobalIndices.

Definition at line 2463 of file Tpetra_CrsGraph_decl.hpp.

template<class LocalOrdinal, class GlobalOrdinal, class Node>
bool Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node >::sortGhostsAssociatedWithEachProcessor_ = true
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 2479 of file Tpetra_CrsGraph_decl.hpp.

Teuchos::RCP<const map_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::map_
protectedinherited

The Map over which this object is distributed.

Definition at line 973 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::imports_
protectedinherited

Buffer into which packed data are imported (received from other processes).

Unfortunately, I had to declare these protected, because CrsMatrix uses them at one point. Please, nobody else use them.

Definition at line 986 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::numImportPacketsPerLID_
protectedinherited

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

Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.

Definition at line 1026 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::exports_
protectedinherited

Buffer from which packed data are exported (sent to other processes).

Unfortunately, I had to declare this protected, because CrsMatrix uses it at one point. Please, nobody else use it.

Definition at line 1033 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< GlobalOrdinal , LocalOrdinal, GlobalOrdinal, Node >::numExportPacketsPerLID_
protectedinherited

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

Unfortunately, I had to declare this protected, because CrsMatrix uses them at one point. Please, nobody else use it.

Definition at line 1048 of file Tpetra_DistObject_decl.hpp.


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