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::Import< LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distribution. More...

#include <Tpetra_Import_decl.hpp>

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

Public Types

using map_type = ::Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node >
 The specialization of Map used by this class. More...
 

Public Member Functions

size_t getNumSameIDs () const
 Number of initial identical IDs. More...
 
size_t getNumPermuteIDs () const
 Number of IDs to permute but not to communicate. More...
 
Kokkos::DualView< const
LocalOrdinal *, device_type > 
getPermuteFromLIDs_dv () const
 List of local IDs in the source Map that are permuted, as a const DualView (that is sync'd to both host and device). More...
 
Teuchos::ArrayView< const
LocalOrdinal > 
getPermuteFromLIDs () const
 List of local IDs in the source Map that are permuted. More...
 
Kokkos::DualView< const
LocalOrdinal *, device_type > 
getPermuteToLIDs_dv () const
 List of local IDs in the target Map that are permuted, as a const DualView (that is sync'd to both host and device). More...
 
Teuchos::ArrayView< const
LocalOrdinal > 
getPermuteToLIDs () const
 List of local IDs in the target Map that are permuted. More...
 
size_t getNumRemoteIDs () const
 Number of entries not on the calling process. More...
 
Kokkos::DualView< const
LocalOrdinal *, device_type > 
getRemoteLIDs_dv () const
 List of entries in the target Map to receive from other processes, as a const DualView (that is sync'd to both host and device). More...
 
Teuchos::ArrayView< const
LocalOrdinal > 
getRemoteLIDs () const
 List of entries in the target Map to receive from other processes. More...
 
size_t getNumExportIDs () const
 Number of entries that must be sent by the calling process to other processes. More...
 
Kokkos::DualView< const
LocalOrdinal *, device_type > 
getExportLIDs_dv () const
 List of entries in the source Map that will be sent to other processes, as a const DualView (that is sync'd to both host and device). More...
 
Teuchos::ArrayView< const
LocalOrdinal > 
getExportLIDs () const
 List of entries in the source Map that will be sent to other processes. More...
 
Teuchos::ArrayView< const int > getExportPIDs () const
 List of processes to which entries will be sent. More...
 
Teuchos::RCP< const map_typegetSourceMap () const
 The source Map used to construct this Export or Import. More...
 
Teuchos::RCP< const map_typegetTargetMap () const
 The target Map used to construct this Export or Import. More...
 
::Tpetra::DistributorgetDistributor () const
 The Distributor that this Export or Import object uses to move data. More...
 
bool isLocallyComplete () const
 Is this Export or Import locally complete? More...
 
bool isLocallyFitted () const
 Are source and target map locally fitted? More...
 
Constructors, assignment, and destructor
 Import (const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target)
 Construct an Import from the source and target Maps. More...
 
 Import (const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target, const Teuchos::RCP< Teuchos::FancyOStream > &out)
 Construct an Import from the source and target Maps, with an output stream for debugging output. More...
 
 Import (const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target, const Teuchos::RCP< Teuchos::ParameterList > &plist)
 Constructor (with list of parameters) More...
 
 Import (const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target, const Teuchos::RCP< Teuchos::FancyOStream > &out, const Teuchos::RCP< Teuchos::ParameterList > &plist)
 Constructor (with list of parameters and debug output stream) More...
 
 Import (const Teuchos::RCP< const map_type > &source, const Teuchos::RCP< const map_type > &target, Teuchos::Array< int > &remotePIDs, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::rcp(new Teuchos::ParameterList))
 Construct an Import from the source and target Maps. More...
 
 Import (const Import< LocalOrdinal, GlobalOrdinal, Node > &importer)
 Copy constructor. More...
 
 Import (const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter)
 "Copy" constructor from an Export object. More...
 
 Import (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &sourceMap, const GlobalOrdinal targetMapRemoteOrPermuteGlobalIndices[], const int targetMapRemoteOrPermuteProcessRanks[], const LocalOrdinal numTargetMapRemoteOrPermuteGlobalIndices, const bool mayReorderTargetMapIndicesLocally, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &out=Teuchos::null)
 Constructor that computes optimized target Map. More...
 
 Import (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &source, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &target, const Teuchos::ArrayView< int > &remotePIDs, const Teuchos::ArrayView< const LocalOrdinal > &userExportLIDs, const Teuchos::ArrayView< const int > &userExportPIDs, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null, const Teuchos::RCP< Teuchos::FancyOStream > &out=Teuchos::null)
 Expert constructor. More...
 
Import< LocalOrdinal,
GlobalOrdinal, Node > & 
operator= (const Import< LocalOrdinal, GlobalOrdinal, Node > &Source)=default
 Assignment operator. More...
 
virtual ~Import ()=default
 Destructor. More...
 
Special "constructors"
void findUnionTargetGIDs (Teuchos::Array< GlobalOrdinal > &unionTgtGIDs, Teuchos::Array< std::pair< int, GlobalOrdinal >> &remotePGIDs, typename Teuchos::Array< GlobalOrdinal >::size_type &numSameGIDs, typename Teuchos::Array< GlobalOrdinal >::size_type &numPermuteGIDs, typename Teuchos::Array< GlobalOrdinal >::size_type &numRemoteGIDs, const Teuchos::ArrayView< const GlobalOrdinal > &sameGIDs1, const Teuchos::ArrayView< const GlobalOrdinal > &sameGIDs2, Teuchos::Array< GlobalOrdinal > &permuteGIDs1, Teuchos::Array< GlobalOrdinal > &permuteGIDs2, Teuchos::Array< GlobalOrdinal > &remoteGIDs1, Teuchos::Array< GlobalOrdinal > &remoteGIDs2, Teuchos::Array< int > &remotePIDs1, Teuchos::Array< int > &remotePIDs2) const
 Find the union of the target IDs from two Import objects. More...
 
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
setUnion (const Import< LocalOrdinal, GlobalOrdinal, Node > &rhs) const
 Return the union of this Import and rhs. More...
 
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
setUnion () const
 Return the union of this Import this->getSourceMap() More...
 
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
createRemoteOnlyImport (const Teuchos::RCP< const map_type > &remoteTarget) const
 Returns an importer that contains only the remote entries of this. More...
 
I/O Methods
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Describe this object in a human-readable way to the given output stream. More...
 
virtual void print (std::ostream &os) const
 Print the Import's data to the given output stream. More...
 

Protected Member Functions

Teuchos::FancyOStream & verboseOutputStream () const
 Valid (nonnull) output stream for verbose output. More...
 
bool verbose () const
 Whether to print verbose debugging output. More...
 
void describeImpl (Teuchos::FancyOStream &out, const std::string &className, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Implementation of describe() for subclasses (Tpetra::Import and Tpetra::Export). More...
 

Protected Attributes

Teuchos::RCP< ImportExportData
< LocalOrdinal, GlobalOrdinal,
Node > > 
TransferData_
 All the data needed for executing the Export communication plan. More...
 

Related Functions

(Note that these are not member functions.)

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
createImport (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &src, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &tgt)
 Nonmember constructor for Import. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
createImport (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &src, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &tgt, const Teuchos::RCP< Teuchos::ParameterList > &plist)
 Nonmember constructor for Import that takes a ParameterList. More...
 

Detailed Description

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

Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distribution.

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.

Tpetra users should use this class to construct a communication plan between two data distributions (i.e., two Map objects). The plan can be called repeatedly by computational classes to perform communication according to the same pattern. Constructing the plan may be expensive, both in terms of communication and computation. However, it can be reused inexpensively.

Tpetra has two classes for data redistribution: Import and Export. Import is for redistributing data from a uniquely-owned distribution to a possibly multiply-owned distribution. Export is for redistributing data from a possibly multiply-owned distribution to a uniquely-owned distribution.

The names "Import" and "Export" have nothing to do with the direction in which data moves relative to the calling process; any process may do both receives and sends in an Import or Export. Rather, the names suggest what happens in their most common use case, the communication pattern for sparse matrix-vector multiply. Import "brings in" remote source vector data (from the domain Map to the column Map) for local computation, and Export "pushes" the result back (from the row Map to the range Map). Import and Export have other uses as well.

As mentioned above, one use case of Import is bringing in remote source vector data for a distributed sparse matrix-vector multiply. The source vector itself is uniquely owned, but must be brought in into an overlapping distribution so that each process can compute its part of the target vector without further communication.

Epetra separated Import and Export for performance reasons. The implementation is different, depending on which direction is the uniquely-owned Map. Tpetra retains this convention.

This class is templated on the same template arguments as Map: the local ordinal type LocalOrdinal, the global ordinal type GlobalOrdinal, and the Kokkos Node type.

This method accepts an optional list of parameters, either through the constructor or through the setParameterList() method. Most users do not need to worry about these parameters; the default values are fine.

Definition at line 109 of file Tpetra_Import_decl.hpp.

Member Typedef Documentation

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

The specialization of Map used by this class.

Definition at line 118 of file Tpetra_Import_decl.hpp.

Constructor & Destructor Documentation

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const map_type > &  source,
const Teuchos::RCP< const map_type > &  target 
)

Construct an Import from the source and target Maps.

Parameters
source[in] The source distribution. This must be a uniquely owned (nonoverlapping) distribution.
target[in] The target distribution. This may be a multiply owned (overlapping) distribution.

Definition at line 188 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const map_type > &  source,
const Teuchos::RCP< const map_type > &  target,
const Teuchos::RCP< Teuchos::FancyOStream > &  out 
)

Construct an Import from the source and target Maps, with an output stream for debugging output.

Parameters
source[in] The source distribution. This must be a uniquely owned (nonoverlapping) distribution.
target[in] The target distribution. This may be a multiply owned (overlapping) distribution.
out[in/out] Output stream for debugging output.

Definition at line 204 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const map_type > &  source,
const Teuchos::RCP< const map_type > &  target,
const Teuchos::RCP< Teuchos::ParameterList > &  plist 
)

Constructor (with list of parameters)

Parameters
source[in] The source distribution. This must be a uniquely owned (nonoverlapping) distribution.
target[in] The target distribution. This may be a multiply owned (overlapping) distribution.
plist[in/out] List of parameters. Currently passed directly to the Distributor that implements communication. If you don't know what this should be, you should use the two-argument constructor, listed above.

Definition at line 215 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const map_type > &  source,
const Teuchos::RCP< const map_type > &  target,
const Teuchos::RCP< Teuchos::FancyOStream > &  out,
const Teuchos::RCP< Teuchos::ParameterList > &  plist 
)

Constructor (with list of parameters and debug output stream)

Parameters
source[in] The source distribution. This must be a uniquely owned (nonoverlapping) distribution.
target[in] The target distribution. This may be a multiply owned (overlapping) distribution.
out[in/out] Output stream (for printing copious debug output on all processes, if that option is enabled).
plist[in/out] List of parameters. Currently passed directly to the Distributor that implements communication. If you don't know what this should be, you should use the two-argument constructor, listed above.

Definition at line 226 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const map_type > &  source,
const Teuchos::RCP< const map_type > &  target,
Teuchos::Array< int > &  remotePIDs,
const Teuchos::RCP< Teuchos::ParameterList > &  plist = Teuchos::rcp(new Teuchos::ParameterList) 
)

Construct an Import from the source and target Maps.

Parameters
source[in] The source distribution. This must be a uniquely owned (nonoverlapping) distribution.
target[in] The target distribution. This may be a multiply owned (overlapping) distribution.
remotePIDs[in] Owning PIDs corresponding to the remoteGIDs. If this information is available one can reduce the cost of the Import constructor.

Definition at line 238 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Import< LocalOrdinal, GlobalOrdinal, Node > &  importer)

Copy constructor.

Note
Currently this only makes a shallow copy of the Import's underlying data.

Definition at line 249 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Export< LocalOrdinal, GlobalOrdinal, Node > &  exporter)

"Copy" constructor from an Export object.

This constructor creates an Import object from the "reverse" of the given Export object. This method is mainly useful for Tpetra developers, for example when building the explicit transpose of a sparse matrix.

Definition at line 255 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  sourceMap,
const GlobalOrdinal  targetMapRemoteOrPermuteGlobalIndices[],
const int  targetMapRemoteOrPermuteProcessRanks[],
const LocalOrdinal  numTargetMapRemoteOrPermuteGlobalIndices,
const bool  mayReorderTargetMapIndicesLocally,
const Teuchos::RCP< Teuchos::ParameterList > &  plist = Teuchos::null,
const Teuchos::RCP< Teuchos::FancyOStream > &  out = Teuchos::null 
)

Constructor that computes optimized target Map.

Like every other Import constructor, this must be called collectively on all processes in the input Map's communicator.

Precondition
sourceMap->isOneToOne()
Parameters
sourceMap[in] Source Map of the Import.
targetMapRemoteOrPermuteGlobalIndices[in] On the calling process, the global indices that will go into the target Map. May differ on each process, just like Map's noncontiguous constructor. No index in here on this process may also appear in sourceMap on this process.
targetMapRemoteOrPermuteProcessRanks[in] For k in 0, ..., numTargetMapRemoteOrPermuteGlobalIndices-1, targetMapRemoteOrPermuteProcessRanks[k] is the rank of the MPI process from which to receive data for global index targetMapRemoteOrPermuteGlobalIndices[k].
numTargetMapRemoteOrPermuteGlobalIndices[in] Number of valid entries in the two input arrays above. May differ on different processes. May be zero on some or even all processes.
mayReorderTargetMapIndicesLocally[in] If true, then this constructor reserves the right to reorder the target Map indices on each process, for better communication performance.

Definition at line 797 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::Import ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  source,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  target,
const Teuchos::ArrayView< int > &  remotePIDs,
const Teuchos::ArrayView< const LocalOrdinal > &  userExportLIDs,
const Teuchos::ArrayView< const int > &  userExportPIDs,
const Teuchos::RCP< Teuchos::ParameterList > &  plist = Teuchos::null,
const Teuchos::RCP< Teuchos::FancyOStream > &  out = Teuchos::null 
)

Expert constructor.

Definition at line 264 of file Tpetra_Import_def.hpp.

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

Destructor.

Member Function Documentation

template<class LocalOrdinal, class GlobalOrdinal, class Node>
Import<LocalOrdinal,GlobalOrdinal,Node>& Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::operator= ( const Import< LocalOrdinal, GlobalOrdinal, Node > &  Source)
default

Assignment operator.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::findUnionTargetGIDs ( Teuchos::Array< GlobalOrdinal > &  unionTgtGIDs,
Teuchos::Array< std::pair< int, GlobalOrdinal >> &  remotePGIDs,
typename Teuchos::Array< GlobalOrdinal >::size_type &  numSameGIDs,
typename Teuchos::Array< GlobalOrdinal >::size_type &  numPermuteGIDs,
typename Teuchos::Array< GlobalOrdinal >::size_type &  numRemoteGIDs,
const Teuchos::ArrayView< const GlobalOrdinal > &  sameGIDs1,
const Teuchos::ArrayView< const GlobalOrdinal > &  sameGIDs2,
Teuchos::Array< GlobalOrdinal > &  permuteGIDs1,
Teuchos::Array< GlobalOrdinal > &  permuteGIDs2,
Teuchos::Array< GlobalOrdinal > &  remoteGIDs1,
Teuchos::Array< GlobalOrdinal > &  remoteGIDs2,
Teuchos::Array< int > &  remotePIDs1,
Teuchos::Array< int > &  remotePIDs2 
) const

Find the union of the target IDs from two Import objects.

On return, the input arrays permuteGIDs[1,2] and remotePGIDs[1,2] will be ordered. unionTgtGIDs are ordered as [{same}, {permute}, {remote}]. {same} is ordered identically to the target map with most "same" indices. {permute} is ordered from smallest ID to largest. {remote} is ordered by remote process ID and ID, respectively. remotePGIDs are ordered the same as {remote}.

Definition at line 1321 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::setUnion ( const Import< LocalOrdinal, GlobalOrdinal, Node > &  rhs) const

Return the union of this Import and rhs.

The "union" of two Import objects is the Import whose source Map is the same as the input Imports' source Maps, and whose target Map is the union of the two input Imports' target Maps. The two input Import objects must have the same source Map. The union operation is symmetric in its two inputs.

The communicator of rhs must be the same as (MPI_EQUAL) or congruent with (MPI_CONGRUENT) this Import's communicator. (That is, the two communicators must have the same number of processes, and each process must have the same rank in both communicators.) This method must be called collectively over that communicator.

The Map that results from this operation does not preserve the original order of global indices in either of the two input Maps. Instead, it sorts global indices on each process so that owned indices occur first, in the same order as in the source Map, and so that remote indices are sorted in order of their sending process rank. This makes communication operations faster.

This primitive is useful for adding two sparse matrices (CrsMatrix), since its can skip over many of the steps of creating the result matrix's column Map from scratch.

We have to call this method "setUnion" rather than "union," because union is a reserved keyword in C++ (and C). It would also be reasonable to call this operator+, though it would be a bit confusing for operator+ to return a pointer but take a reference (or to take a pointer, but have the left-hand side of the + expression be a reference).

Definition at line 1426 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::setUnion ( ) const

Return the union of this Import this->getSourceMap()

This special case of setUnion creates a new Import object such that the targetMap of the new object contains all local unknowns in the sourceMap (plus whatever remotes were contained in this->getSourceMap()).

The Map that results from this operation does not preserve the input order of global indices. All local global indices are ordered in the order of the sourceMap, all remotes are ordered as implied by the Importer for *this.

This primitive is useful for adding or multipyling two sparse matrices (CrsMatrix), since its can skip over many of the steps of creating the result matrix's column Map from scratch.

Definition at line 1683 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::createRemoteOnlyImport ( const Teuchos::RCP< const map_type > &  remoteTarget) const

Returns an importer that contains only the remote entries of this.

Returns an importer that contains only the remote entries of this importer. It is expected that remoteTarget represents such a map.

Definition at line 1754 of file Tpetra_Import_def.hpp.

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

Describe this object in a human-readable way to the given output stream.

You must call this method as a collective over all processes in the communicator of the source and target Map of this object.

Parameters
out[out] Output stream to which to write. Only Process 0 in this object's communicator may write to the output stream.
verbLevel[in] Verbosity level. This also controls whether this method does any communication. At verbosity levels higher (greater) than Teuchos::VERB_LOW, this method behaves as a collective over the object's communicator.

Teuchos::FancyOStream wraps std::ostream. It adds features like tab levels. If you just want to wrap std::cout, try this:

auto out = Teuchos::getFancyOStream (Teuchos::rcpFromRef (std::out));

Reimplemented from Tpetra::Details::Transfer< LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 938 of file Tpetra_Import_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
void Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node >::print ( std::ostream &  os) const
virtual

Print the Import's data to the given output stream.

This method assumes that the given output stream can be written on all process(es) in the Import's communicator. The resulting output is useful mainly for debugging.

Note
This method tries its best (by using barriers at the end of each iteration of a for loop over all communicator ranks) to ensure ordered deterministic output. However, the assumption that all processes can write to the stream means that there are no ordering guarantees other than what the operating and run-time system provide. (MPI synchronization may be separate from output stream synchronization, so the barriers only improve the chances that output can complete before the next process starts writing.)

Definition at line 947 of file Tpetra_Import_def.hpp.

size_t Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getNumSameIDs ( ) const
inherited

Number of initial identical IDs.

The number of IDs that are identical between the source and target Maps, up to the first different ID.

size_t Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getNumPermuteIDs ( ) const
inherited

Number of IDs to permute but not to communicate.

The number of IDs that are local to the calling process, but not part of the first getNumSameIDs() entries. The Import or Export will permute these entries locally (without distributed-memory communication).

Kokkos::DualView<const LocalOrdinal *, device_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getPermuteFromLIDs_dv ( ) const
inherited

List of local IDs in the source Map that are permuted, as a const DualView (that is sync'd to both host and device).

Teuchos::ArrayView<const LocalOrdinal > Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getPermuteFromLIDs ( ) const
inherited

List of local IDs in the source Map that are permuted.

Kokkos::DualView<const LocalOrdinal *, device_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getPermuteToLIDs_dv ( ) const
inherited

List of local IDs in the target Map that are permuted, as a const DualView (that is sync'd to both host and device).

Teuchos::ArrayView<const LocalOrdinal > Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getPermuteToLIDs ( ) const
inherited

List of local IDs in the target Map that are permuted.

size_t Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getNumRemoteIDs ( ) const
inherited

Number of entries not on the calling process.

Kokkos::DualView<const LocalOrdinal *, device_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getRemoteLIDs_dv ( ) const
inherited

List of entries in the target Map to receive from other processes, as a const DualView (that is sync'd to both host and device).

Teuchos::ArrayView<const LocalOrdinal > Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getRemoteLIDs ( ) const
inherited

List of entries in the target Map to receive from other processes.

size_t Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getNumExportIDs ( ) const
inherited

Number of entries that must be sent by the calling process to other processes.

Kokkos::DualView<const LocalOrdinal *, device_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getExportLIDs_dv ( ) const
inherited

List of entries in the source Map that will be sent to other processes, as a const DualView (that is sync'd to both host and device).

Teuchos::ArrayView<const LocalOrdinal > Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getExportLIDs ( ) const
inherited

List of entries in the source Map that will be sent to other processes.

Teuchos::ArrayView<const int> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getExportPIDs ( ) const
inherited

List of processes to which entries will be sent.

The entry with local ID getExportLIDs()[i] will be sent to process getExportPiDs()[i].

Teuchos::RCP<const map_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getSourceMap ( ) const
inherited

The source Map used to construct this Export or Import.

Teuchos::RCP<const map_type> Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getTargetMap ( ) const
inherited

The target Map used to construct this Export or Import.

::Tpetra::Distributor& Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::getDistributor ( ) const
inherited

The Distributor that this Export or Import object uses to move data.

bool Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::isLocallyComplete ( ) const
inherited

Is this Export or Import locally complete?

If this is an Export, then do all source Map indices on the calling process exist on at least one process (not necessarily this one) in the target Map?

If this is an Import, then do all target Map indices on the calling process exist on at least one process (not necessarily this one) in the source Map?

It's not necessarily an error for an Export or Import not to be locally complete on one or more processes. For example, this may happen in the common use case of "restriction" – that is, taking a subset of a large object. Nevertheless, you may find this predicate useful for figuring out whether you set up your Maps in the way that you expect.

bool Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::isLocallyFitted ( ) const
inherited

Are source and target map locally fitted?

Returns whether source and target map are locally fitted on the calling rank. This is can be more efficient that calling isLocallyFitted() on the maps directly, since no indices need to be compared.

Teuchos::FancyOStream& Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::verboseOutputStream ( ) const
protectedinherited

Valid (nonnull) output stream for verbose output.

bool Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::verbose ( ) const
protectedinherited

Whether to print verbose debugging output.

void Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::describeImpl ( Teuchos::FancyOStream &  out,
const std::string &  className,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
protectedinherited

Implementation of describe() for subclasses (Tpetra::Import and Tpetra::Export).

Parameters
out[out] Output stream to which to write. Only Process 0 in this object's communicator may write to the output stream.
className[in] Name of the subclass of Transfer calling this method.
verbLevel[in] Verbosity level. This also controls whether this method does any communication. At verbosity levels higher (greater) than Teuchos::VERB_LOW, this method behaves as a collective over the object's communicator.

Friends And Related Function Documentation

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > createImport ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  src,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  tgt 
)
related

Nonmember constructor for Import.

Create a Import object from the given source and target Maps.

Precondition
src != null
tgt != null
Returns
The Import object. If src == tgt, returns null. (Debug mode: throws std::runtime_error if one of src or tgt is null.)

Definition at line 532 of file Tpetra_Import_decl.hpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > createImport ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  src,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  tgt,
const Teuchos::RCP< Teuchos::ParameterList > &  plist 
)
related

Nonmember constructor for Import that takes a ParameterList.

Create a Import object from the given source and target Maps, using the given list of parameters.

Precondition
src != null
tgt != null
Returns
The Import object. If src == tgt, returns null. (Debug mode: throws std::runtime_error if one of src or tgt is null.)

Definition at line 559 of file Tpetra_Import_decl.hpp.

Member Data Documentation

Teuchos::RCP<ImportExportData<LocalOrdinal , GlobalOrdinal , Node > > Tpetra::Details::Transfer< LocalOrdinal , GlobalOrdinal , Node >::TransferData_
protectedinherited

All the data needed for executing the Export communication plan.

Definition at line 271 of file Tpetra_Details_Transfer_decl.hpp.


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