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 Types | Protected Member Functions | List of all members
Tpetra::Details::CooMatrix< SC, LO, GO, NT > Class Template Referenceabstract

Sparse matrix used only for file input / output. More...

#include <Tpetra_Details_CooMatrix.hpp>

Inheritance diagram for Tpetra::Details::CooMatrix< SC, LO, GO, NT >:
Inheritance graph
[legend]

Public Types

typedef char packet_type
 This class transfers data as bytes (MPI_BYTE). More...
 
typedef SC scalar_type
 Type of each entry (value) in the sparse matrix. More...
 
typedef ::Tpetra::Map
< local_ordinal_type,
global_ordinal_type, node_type > 
map_type
 Type of the Map specialization to give to the constructor. More...
 
Typedefs
using execution_space = typename device_type::execution_space
 The Kokkos execution space. More...
 

Public Member Functions

 CooMatrix ()
 Default constructor. More...
 
 CooMatrix (const ::Teuchos::RCP< const map_type > &map)
 Constructor that takes a Map. More...
 
virtual ~CooMatrix ()
 Destructor (virtual for memory safety of derived classes). More...
 
void sumIntoGlobalValue (const GO gblRowInd, const GO gblColInd, const SC &val)
 Insert one entry locally into the sparse matrix, if it does not exist there yet. If it does exist, sum the values. More...
 
void sumIntoGlobalValues (const GO gblRowInds[], const GO gblColInds[], const SC vals[], const std::size_t numEnt)
 Insert multiple entries locally into the sparse matrix. More...
 
void sumIntoGlobalValues (std::initializer_list< GO > gblRowInds, std::initializer_list< GO > gblColInds, std::initializer_list< SC > vals, const std::size_t numEnt)
 Initializer-list overload of the above method (which see). More...
 
std::size_t getLclNumEntries () const
 Number of entries in the sparse matrix on the calling process. More...
 
void fillComplete (const ::Teuchos::RCP< const ::Teuchos::Comm< int > > &comm)
 Tell the matrix that you are done inserting entries locally, and that the matrix should build its Map now. More...
 
void fillComplete ()
 Special version of fillComplete that assumes that the matrix already has a Map, and reuses its communicator to create a new Map. More...
 
bool localError () const
 Whether this object had an error on the calling process. More...
 
std::string errorMessages () const
 The current stream of error messages. More...
 
virtual std::string description () const
 One-line descriptiion of this object; overrides Teuchos::Describable method. More...
 
virtual void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print a descriptiion of this object to the given output stream; overrides Teuchos::Describable method. More...
 
Public methods for redistributing data
void doImport (const SrcDistObject &source, const Import< LO, GO, NT > &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< LO, GO, NT > &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< LO, GO, NT > &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< LO, GO, NT > &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< LO, GO, NT > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void beginImport (const SrcDistObject &source, const Export< LO, GO, NT > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Export< LO, GO, NT > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void beginExport (const SrcDistObject &source, const Import< LO, GO, NT > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Import< LO, GO, NT > &importer, const CombineMode CM, const bool restrictedMode=false)
 
void endImport (const SrcDistObject &source, const Export< LO, GO, NT > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Export< LO, GO, NT > &exporter, const CombineMode CM, const bool restrictedMode=false)
 
void endExport (const SrcDistObject &source, const Import< LO, GO, NT > &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...
 

Protected Types

using buffer_device_type = typename::Tpetra::DistObject< char, LO, GO, NT >::buffer_device_type
 Kokkos::Device specialization for DistObject communication buffers. More...
 

Protected Member Functions

virtual size_t constantNumberOfPackets () const
 By returning 0, tell DistObject that this class may not necessarily have a constant number of "packets" per local index. More...
 
virtual bool checkSizes (const ::Tpetra::SrcDistObject &source)
 Compare the source and target (this) objects for compatibility. 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 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.")

virtual bool checkSizes (const SrcDistObject &source)=0
 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)
 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...
 
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...
 

Detailed Description

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
class Tpetra::Details::CooMatrix< SC, LO, GO, NT >

Sparse matrix used only for file input / output.

This class stores a sparse matrix in coordinate format. It is meant only to help file input and output. Thus, it implements DistObject, but does NOT implement RowMatrix or even Operator.

Unlike other DistObject subclasses in Tpetra, this class' constructor need not necessarily take a Map. If the class does NOT have a Map, it builds its Map at fillComplete(), after construction, as a function of the input indices.

Users are only allowed to insert matrix entries if the class does NOT have a Map. Insertion is local to each process. Users insert entries by calling sumIntoGlobalValue(). Each entry is a "triple" consisting of a global row index, a global column index, and a matrix value. Users call fillComplete() when they are done inserting entries. At that point, the object builds its Map. (Unlike CrsMatrix, this class does NOT allow multiple resumeFill() / fillComplete() cycles.)

Once the class has a Map, users may apply DistObject methods like doImport() and doExport() to redistribute the data. The target of an Import or Export must have a Map, and users may not have inserted entries into it.

Here is an example of how to use this class:

for (size_t k = 0; k < numEntriesToInsert; ++k) {
long long gblRowInd, gblColInd;
double val;
// We don't implement this here. You have to do it.
readEntryFromFile (&gblRowInd, &gblColInd, &val);
A_in.sumIntoGlobalValue (gblRowInd, gblColInd, val);
}
// You are responsible for supplying the output matrix's
// communicator and Map.
A_in.fillComplete (comm);
Teuchos::RCP<const Tpetra::Map<int, long long> > outMap = ...;
Tpetra::Export<int, long long> exporter (A_in.getMap (), outMap);
A_out.doExport (A_in, exporter, Tpetra::ADD);

Definition at line 574 of file Tpetra_Details_CooMatrix.hpp.

Member Typedef Documentation

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
typedef char Tpetra::Details::CooMatrix< SC, LO, GO, NT >::packet_type

This class transfers data as bytes (MPI_BYTE).

Definition at line 577 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
typedef SC Tpetra::Details::CooMatrix< SC, LO, GO, NT >::scalar_type

Type of each entry (value) in the sparse matrix.

Definition at line 579 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
typedef ::Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> Tpetra::Details::CooMatrix< SC, LO, GO, NT >::map_type

Type of the Map specialization to give to the constructor.

Definition at line 585 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
using Tpetra::Details::CooMatrix< SC, LO, GO, NT >::buffer_device_type = typename ::Tpetra::DistObject<char, LO, GO, NT>::buffer_device_type
protected

Kokkos::Device specialization for DistObject communication buffers.

Definition at line 999 of file Tpetra_Details_CooMatrix.hpp.

using Tpetra::DistObject< char , LO , GO , NT >::execution_space = typename device_type::execution_space
inherited

The Kokkos execution space.

Definition at line 345 of file Tpetra_DistObject_decl.hpp.

Constructor & Destructor Documentation

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
Tpetra::Details::CooMatrix< SC, LO, GO, NT >::CooMatrix ( )
inline

Default constructor.

This creates the object with a null Map. Users may insert entries into this object, until they call fillComplete. This object may NOT be the target of an Import or Export operation.

Definition at line 601 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
Tpetra::Details::CooMatrix< SC, LO, GO, NT >::CooMatrix ( const ::Teuchos::RCP< const map_type > &  map)
inline

Constructor that takes a Map.

Parameters
map[in] Input Map.

If map is nonnull, then users may use this object ONLY as the target of an Import or Export operation. They may NOT insert entries into it.

Definition at line 614 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
virtual Tpetra::Details::CooMatrix< SC, LO, GO, NT >::~CooMatrix ( )
inlinevirtual

Destructor (virtual for memory safety of derived classes).

Definition at line 621 of file Tpetra_Details_CooMatrix.hpp.

Member Function Documentation

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::sumIntoGlobalValue ( const GO  gblRowInd,
const GO  gblColInd,
const SC &  val 
)
inline

Insert one entry locally into the sparse matrix, if it does not exist there yet. If it does exist, sum the values.

Parameters
gblRowInd[in] Global row index of the entry to insert.
gblColInd[in] Global column index of the entry to insert.
val[in] Value of the matrix entry to insert / sum.

Definition at line 630 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::sumIntoGlobalValues ( const GO  gblRowInds[],
const GO  gblColInds[],
const SC  vals[],
const std::size_t  numEnt 
)
inline

Insert multiple entries locally into the sparse matrix.

This works like multiple calls to sumIntoGlobalValue.

Parameters
gblRowInd[in] Global row indices of the entries to insert.
gblColInd[in] Global column indices of the entries to insert.
val[in] Values of the matrix entries to insert / sum.
numEnt[in] Number of entries to insert.

Definition at line 646 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::sumIntoGlobalValues ( std::initializer_list< GO >  gblRowInds,
std::initializer_list< GO >  gblColInds,
std::initializer_list< SC >  vals,
const std::size_t  numEnt 
)
inline

Initializer-list overload of the above method (which see).

Definition at line 656 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
std::size_t Tpetra::Details::CooMatrix< SC, LO, GO, NT >::getLclNumEntries ( ) const
inline

Number of entries in the sparse matrix on the calling process.

Definition at line 667 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::fillComplete ( const ::Teuchos::RCP< const ::Teuchos::Comm< int > > &  comm)
inline

Tell the matrix that you are done inserting entries locally, and that the matrix should build its Map now.

This is the preferred version of fillComplete().

Precondition
The matrix does not yet have a Map.
fillComplete() has never been called before on this object.
Postcondition
The matrix has a Map.
Parameters
comm[in] Input communicator to use for the Map.

Definition at line 717 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::fillComplete ( )
inline

Special version of fillComplete that assumes that the matrix already has a Map, and reuses its communicator to create a new Map.

DO NOT call this method unless you know what you are doing.

Precondition
The matrix DOES have a (nonnull) Map with a nonnull communicator.

Definition at line 736 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
bool Tpetra::Details::CooMatrix< SC, LO, GO, NT >::localError ( ) const
inline

Whether this object had an error on the calling process.

Import and Export operations using this object as the target of the Import or Export may incur local errors, for example due to invalid rows or incorrectly sized data. On local error detection, we don't want to throw an exception right away, because not all processes may throw an exception; this inconsistency across processes can result in deadlock or put Tpetra in an incorrect state. Instead, we set a local error flag on the affected process, and ignore the incorrect data. If you want to check whether any process experienced an error, you must do a reduction or all-reduce over this flag. Every time you initiate a new Import or Export with this object as the target, we clear this flag.

Definition at line 760 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
std::string Tpetra::Details::CooMatrix< SC, LO, GO, NT >::errorMessages ( ) const
inline

The current stream of error messages.

This is only nonempty on the calling process if localError() returns true. In that case, it stores a stream of human-readable, endline-separated error messages encountered during an Import or Export cycle. Every time you initiate a new Import or Export with this object as the target, we clear this stream.

If you want to print this, you are responsible for ensuring that it is valid for the calling MPI process to print to whatever output stream you use. On some MPI implementations, you may need to send the string to Process 0 in MPI_COMM_WORLD for printing.

Note to developers: we clear the stream at the beginning of checkSizes().

Definition at line 781 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
virtual std::string Tpetra::Details::CooMatrix< SC, LO, GO, NT >::description ( ) const
inlinevirtual

One-line descriptiion of this object; overrides Teuchos::Describable method.

Reimplemented from Tpetra::DistObject< char, LO, GO, NT >.

Definition at line 823 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
virtual void Tpetra::Details::CooMatrix< SC, LO, GO, NT >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
inlinevirtual

Print a descriptiion of this object to the given output stream; overrides Teuchos::Describable method.

Reimplemented from Tpetra::DistObject< char, LO, GO, NT >.

Definition at line 843 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
virtual size_t Tpetra::Details::CooMatrix< SC, LO, GO, NT >::constantNumberOfPackets ( ) const
inlineprotectedvirtual

By returning 0, tell DistObject that this class may not necessarily have a constant number of "packets" per local index.

Reimplemented from Tpetra::DistObject< char, LO, GO, NT >.

Definition at line 967 of file Tpetra_Details_CooMatrix.hpp.

template<class SC, class LO = ::Tpetra::DistObject<char>::local_ordinal_type, class GO = ::Tpetra::DistObject<char>::global_ordinal_type, class NT = ::Tpetra::DistObject<char>::node_type>
virtual bool Tpetra::Details::CooMatrix< SC, LO, GO, NT >::checkSizes ( const ::Tpetra::SrcDistObject source)
inlineprotectedvirtual

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

Returns
True if they are compatible, else false.

Definition at line 975 of file Tpetra_Details_CooMatrix.hpp.

void Tpetra::DistObject< char , LO , GO , NT >::doImport ( const SrcDistObject< char, LO, GO, NT > &  source,
const Import< LO , GO , NT > &  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< char , LO , GO , NT >::doImport ( const SrcDistObject< char, LO, GO, NT > &  source,
const Export< LO , GO , NT > &  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< char , LO , GO , NT >::doExport ( const SrcDistObject< char, LO, GO, NT > &  source,
const Export< LO , GO , NT > &  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< char , LO , GO , NT >::doExport ( const SrcDistObject< char, LO, GO, NT > &  source,
const Import< LO , GO , NT > &  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< char , LO , GO , NT >::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< char , LO , GO , NT >::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< char , LO , GO , NT >::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 589 of file Tpetra_DistObject_decl.hpp.

void Tpetra::DistObject< char , LO , GO , NT >::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< char , LO , GO , NT >::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.
virtual void Tpetra::DistObject< char , LO , GO , NT >::doTransfer ( const SrcDistObject< char, LO, GO, NT > &  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< char , LO , GO , NT >::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< char , LO , GO , NT >::beginTransfer ( const SrcDistObject< char, LO, GO, NT > &  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 bool Tpetra::DistObject< char , LO , GO , NT >::checkSizes ( const SrcDistObject< char, LO, GO, NT > &  source)
protectedpure virtualinherited

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

Returns
True if they are compatible, else false.
virtual void Tpetra::DistObject< char , LO , GO , NT >::copyAndPermute ( const SrcDistObject< char, LO, GO, NT > &  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< char , LO , GO , NT >::copyAndPermute ( const SrcDistObject< char, LO, GO, NT > &  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< char , LO , GO , NT >::packAndPrepare ( const SrcDistObject< char, LO, GO, NT > &  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< char , LO , GO , NT >::packAndPrepare ( const SrcDistObject< char, LO, GO, NT > &  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< char , LO , GO , NT >::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< char , LO , GO , NT >::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.

Member Data Documentation

Teuchos::RCP<const map_type> Tpetra::DistObject< char , LO , GO , NT >::map_
protectedinherited

The Map over which this object is distributed.

Definition at line 998 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< char , LO , GO , NT >::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 1011 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< char , LO , GO , NT >::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 1051 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<packet_type*, buffer_device_type> Tpetra::DistObject< char , LO , GO , NT >::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 1058 of file Tpetra_DistObject_decl.hpp.

Kokkos::DualView<size_t*, buffer_device_type> Tpetra::DistObject< char , LO , GO , NT >::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 1073 of file Tpetra_DistObject_decl.hpp.


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