Tpetra parallel linear algebra  Version of the Day
Namespaces | Classes | Enumerations | Functions
Tpetra::Details Namespace Reference

Namespace for Tpetra implementation details. More...

Namespaces

 DefaultTypes
 Declarations of values of Tpetra classes' default template parameters.
 

Classes

struct  AbsMax
 Functor for the the ABSMAX CombineMode of Import and Export operations. More...
 
class  CommRequest
 Base class for the request (more or less a future) representing a pending nonblocking MPI operation. More...
 
class  ContiguousUniformDirectory
 Implementation of Directory for a contiguous, uniformly distributed Map. More...
 
class  CooMatrix
 Sparse matrix used only for file input / output. More...
 
struct  CrsMatrixGetDiagCopyFunctor
 Functor that implements much of the one-argument overload of Tpetra::CrsMatrix::getLocalDiagCopy, for the case where the matrix is fill complete. More...
 
class  Directory
 Computes the local ID and process ID corresponding to given global IDs. More...
 
class  DistributedContiguousDirectory
 Implementation of Directory for a distributed contiguous Map. More...
 
class  DistributedNoncontiguousDirectory
 Implementation of Directory for a distributed noncontiguous Map. More...
 
class  FixedHashTable
 
struct  Hash
 The hash function for FixedHashTable. More...
 
struct  Hash< KeyType, DeviceType, OffsetType, int >
 Specialization for ResultType = int. More...
 
class  HashTable
 
class  InvalidGlobalIndex
 Exception thrown by CrsMatrix on invalid global index. More...
 
class  InvalidGlobalRowIndex
 Exception thrown by CrsMatrix on invalid global row index. More...
 
class  LocalMap
 "Local" part of Map suitable for Kokkos kernels. More...
 
class  MapCloner
 Implementation detail of Map::clone(). More...
 
struct  MultiVectorCloner
 Implementation of Tpetra::MultiVector::clone(). More...
 
class  MultiVectorFillerData
 Implementation of fill and local assembly for MultiVectorFiller. More...
 
class  MultiVectorFillerData2
 Second implementation of fill and local assembly for MultiVectorFiller. More...
 
class  OptColMap
 Implementation detail of makeOptimizedColMap, and makeOptimizedColMapAndImport. More...
 
struct  PackCrsMatrixError
 Reduction result for PackCrsMatrixFunctor below. More...
 
struct  PackTraits
 Traits class for packing / unpacking data of type T, using Kokkos data structures that live in the given space D. More...
 
class  ProfilingRegion
 Profile the given scope. More...
 
class  ReplicatedDirectory
 Implementation of Directory for a locally replicated Map. More...
 
class  TieBreak
 Interface for breaking ties in ownership. More...
 
class  Transfer
 Common base class of Import and Export. More...
 
struct  UnpackCrsMatrixAndCombineFunctor
 Unpacks and combines a single row of the CrsMatrix. More...
 
struct  UnpackCrsMatrixError
 Reduction result for UnpackCrsMatrixAndCombineFunctor below. More...
 

Enumerations

enum  EStorageStatus
 Status of the graph's or matrix's storage, when not in a fill-complete state. More...
 
enum  EDistributorSendType
 The type of MPI send that Distributor should use. More...
 
enum  EDistributorHowInitialized
 Enum indicating how and whether a Distributor was initialized. More...
 

Functions

template<class LO , class GO , class DT , class OffsetType , class NumEntType >
OffsetType convertColumnIndicesFromGlobalToLocal (const Kokkos::View< LO *, DT > &lclColInds, const Kokkos::View< const GO *, DT > &gblColInds, const Kokkos::View< const OffsetType *, DT > &ptr, const LocalMap< LO, GO, DT > &lclColMap, const Kokkos::View< const NumEntType *, DT > &numRowEnt)
 Convert a (StaticProfile) CrsGraph's global column indices into local column indices. More...
 
template<class OffsetsViewType , class CountsViewType , class SizeType = typename OffsetsViewType::size_type>
OffsetsViewType::non_const_value_type computeOffsetsFromCounts (const OffsetsViewType &ptr, const CountsViewType &counts)
 Compute offsets from counts. More...
 
template<class OffsetsViewType , class CountType , class SizeType = typename OffsetsViewType::size_type>
OffsetsViewType::non_const_value_type computeOffsetsFromConstantCount (const OffsetsViewType &ptr, const CountType &count)
 Compute offsets from a constant count. More...
 
template<class OutputViewType , class InputViewType >
void copyConvert (const OutputViewType &dst, const InputViewType &src)
 Copy values from the 1-D Kokkos::View src, to the 1-D Kokkos::View dst, of the same length. The entries of src and dst may have different types, but it must be possible to copy-construct each entry of dst with its corresponding entry of src. More...
 
template<class OutputViewType , class InputViewType >
void copyOffsets (const OutputViewType &dst, const InputViewType &src)
 Copy row offsets (in a sparse graph or matrix) from src to dst. The offsets may have different types. More...
 
template<class ValueType , class OutputDeviceType >
Impl::CreateMirrorViewFromUnmanagedHostArray< ValueType, OutputDeviceType >::output_view_type create_mirror_view_from_raw_host_array (const OutputDeviceType &, ValueType *inPtr, const size_t inSize, const bool copy=true, const char label[]="")
 Variant of Kokkos::create_mirror_view that takes a raw host 1-d array as input. More...
 
template<class SparseMatrixType , class ValsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type crsMatrixSumIntoValues_sortedSortedLinear (const SparseMatrixType &A, const typename SparseMatrixType::ordinal_type lclRow, const typename SparseMatrixType::ordinal_type lclColInds[], const typename SparseMatrixType::ordinal_type sortPerm[], const ValsViewType &vals, const typename SparseMatrixType::ordinal_type numEntInInput, const bool forceAtomic=false, const bool checkInputIndices=true)
  A(lclRow, lclColsInds[sortPerm[j]]) += vals[sortPerm[j]], for all j in 0 .. eltDim-1. More...
 
template<class SparseMatrixType , class ValsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type crsMatrixReplaceValues_sortedSortedLinear (const SparseMatrixType &A, const typename SparseMatrixType::ordinal_type lclRow, const typename SparseMatrixType::ordinal_type lclColInds[], const typename SparseMatrixType::ordinal_type sortPerm[], const ValsViewType &vals, const typename SparseMatrixType::ordinal_type numEntInInput, const bool forceAtomic=false, const bool checkInputIndices=true)
  A(lclRow, lclColsInds[sortPerm[j]]) = vals[sortPerm[j]], for all j in 0 .. eltDim-1. More...
 
template<class SparseMatrixType , class VectorViewType , class RhsViewType , class LhsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type crsMatrixAssembleElement_sortedLinear (const SparseMatrixType &A, const VectorViewType &x, typename SparseMatrixType::ordinal_type lids[], typename SparseMatrixType::ordinal_type sortPerm[], const RhsViewType &rhs, const LhsViewType &lhs, const bool forceAtomic=false, const bool checkInputIndices=true)
 A(lids[j], lids[j]) += lhs(j,j) and x(lids[j]) += rhs(j), for all j in 0 .. eltDim-1. More...
 
void gathervPrint (std::ostream &out, const std::string &s, const Teuchos::Comm< int > &comm)
 On Process 0 in the given communicator, print strings from each process in that communicator, in rank order. More...
 
template<class DiagType , class LocalMapType , class CrsMatrixType >
static LocalMapType::local_ordinal_type getDiagCopyWithoutOffsets (const DiagType &D, const LocalMapType &rowMap, const LocalMapType &colMap, const CrsMatrixType &A)
 Given a locally indexed, local sparse matrix, and corresponding local row and column Maps, extract the matrix's diagonal entries into a 1-D Kokkos::View. More...
 
template<class SC , class LO , class GO , class NT >
LO getLocalDiagCopyWithoutOffsetsNotFillComplete (::Tpetra::Vector< SC, LO, GO, NT > &diag, const ::Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool debug=false)
 Given a locally indexed, global sparse matrix, extract the matrix's diagonal entries into a Tpetra::Vector. More...
 
template<class InputViewType , class OutputViewType >
std::shared_ptr< CommRequestiallreduce (const InputViewType &sendbuf, const OutputViewType &recvbuf, const ::Teuchos::EReductionType op, const ::Teuchos::Comm< int > &comm)
 Nonblocking all-reduce, for either rank-1 or rank-0 Kokkos::View objects. More...
 
bool isInterComm (const Teuchos::Comm< int > &comm)
 Return true if and only if the input communicator wraps an MPI intercommunicator. More...
 
template<class LO , class GO , class NT >
int makeColMap (Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &colMap, const Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &domMap, const RowGraph< LO, GO, NT > &graph, const bool sortEachProcsGids=true, std::ostream *errStrm=NULL)
 Make the graph's column Map. More...
 
template<class MapType >
MapType makeOptimizedColMap (std::ostream &errStream, bool &lclErr, const MapType &domMap, const MapType &colMap)
 Return an optimized reordering of the given column Map. More...
 
template<class MapType >
std::pair< MapType, Teuchos::RCP< typename OptColMap< MapType >::import_type > > makeOptimizedColMapAndImport (std::ostream &errStream, bool &lclErr, const MapType &domMap, const MapType &colMap, const typename OptColMap< MapType >::import_type *oldImport, const bool makeImport)
 Return an optimized reordering of the given column Map. Optionally, recompute an Import from the input domain Map to the new column Map. More...
 
template<class OrdinalType , class IndexType >
IndexType countMergeUnsortedIndices (const OrdinalType curInds[], const IndexType numCurInds, const OrdinalType inputInds[], const IndexType numInputInds)
 Count the number of column indices that can be merged into the current row, assuming that both the current row's indices and the input indices are unsorted. More...
 
template<class OrdinalType , class IndexType >
IndexType countMergeSortedIndices (const OrdinalType curInds[], const IndexType numCurInds, const OrdinalType inputInds[], const IndexType numInputInds)
 Count the number of column indices that can be merged into the current row, assuming that both the current row's indices and the input indices are sorted. More...
 
template<class OrdinalType , class IndexType >
std::pair< bool, IndexType > mergeSortedIndices (OrdinalType curInds[], const IndexType midPos, const IndexType endPos, const OrdinalType inputInds[], const IndexType numInputInds)
 Attempt to merge the input indices into the current row's column indices, assuming that both the current row's indices and the input indices are sorted. More...
 
template<class OrdinalType , class IndexType >
std::pair< bool, IndexType > mergeUnsortedIndices (OrdinalType curInds[], const IndexType midPos, const IndexType endPos, const OrdinalType inputInds[], const IndexType numInputInds)
 Attempt to merge the input indices into the current row's column indices, assuming that both the current row's indices and the input indices are unsorted. More...
 
template<class OrdinalType , class ValueType , class IndexType >
std::pair< bool, IndexType > mergeUnsortedIndicesAndValues (OrdinalType curInds[], ValueType curVals[], const IndexType midPos, const IndexType endPos, const OrdinalType inputInds[], const ValueType inputVals[], const IndexType numInputInds)
 Attempt to merge the input indices and values into the current row's column indices and corresponding values, assuming that both the current row's indices and the input indices are unsorted. More...
 
template<class LocalMatrixType , class LocalMapType >
KOKKOS_FUNCTION bool packCrsMatrixRow (const LocalMatrixType &lclMatrix, const LocalMapType &lclColMap, char *const numEntOut, char *const valOut, char *const indOut, const size_t numEnt, const typename LocalMatrixType::ordinal_type lclRow)
 Packs a single row of the CrsMatrix. More...
 
template<class LocalMatrixType , class LocalMapType >
bool packCrsMatrix (const LocalMatrixType &lclMatrix, const LocalMapType &lclColMap, std::unique_ptr< std::string > &errStr, Teuchos::Array< char > &exports, const Teuchos::ArrayView< size_t > &numPacketsPerLID, size_t &constantNumPackets, const Teuchos::ArrayView< const typename LocalMatrixType::ordinal_type > &exportLIDs, const int myRank, Distributor &)
 Pack specified entries of the given local sparse matrix for communication. More...
 
template<class KeyType , class ValueType >
KOKKOS_FUNCTION void shortSortKeysAndValues_2 (KeyType keys[2], ValueType values[2])
 Sort keys and values jointly, by keys, for arrays of length 2. More...
 
template<class KeyType >
KOKKOS_FUNCTION void shortSortKeys_2 (KeyType keys[2])
 Sort length-2 array of keys. More...
 
template<class KeyType , class ValueType >
KOKKOS_FUNCTION void shortSortKeysAndValues_3 (KeyType keys[3], ValueType values[3])
 Sort keys and values jointly, by keys, for arrays of length 3. More...
 
template<class KeyType >
KOKKOS_FUNCTION void shortSortKeys_3 (KeyType keys[3])
 Sort length-3 array of keys. More...
 
template<class KeyType , class ValueType >
KOKKOS_FUNCTION void shortSortKeysAndValues_4 (KeyType keys[4], ValueType values[4])
 Sort keys and values jointly, by keys, for arrays of length 4. More...
 
template<class KeyType >
KOKKOS_FUNCTION void shortSortKeys_4 (KeyType keys[4])
 Sort length-4 array of keys. More...
 
template<class KeyType , class ValueType >
KOKKOS_FUNCTION void shortSortKeysAndValues_8 (KeyType keys[8], ValueType values[8])
 Sort keys and values jointly, by keys, for arrays of length 8. More...
 
template<class KeyType >
KOKKOS_FUNCTION void shortSortKeys_8 (KeyType keys[8])
 Sort length-8 array of keys. More...
 
template<class KeyType , class ValueType , class IndexType >
KOKKOS_FUNCTION void shellSortKeysAndValues (KeyType keys[], ValueType values[], const IndexType n)
 Shellsort (yes, it's one word) the input array keys, and apply the resulting permutation to the input array values. More...
 
template<class KeyType , class IndexType >
KOKKOS_FUNCTION void shellSortKeys (KeyType keys[], const IndexType n)
 Shellsort (yes, it's one word) the input array keys. More...
 
template<class LocalMatrixType , class LocalMapType >
bool unpackCrsMatrixAndCombine (LocalMatrixType &lclMatrix, const LocalMapType &lclColMap, std::unique_ptr< std::string > &errStr, const Teuchos::ArrayView< const typename LocalMatrixType::ordinal_type > &importLIDs, const Teuchos::ArrayView< const char > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, size_t constantNumPackets, const int myRank, Distributor &, CombineMode combineMode, const bool atomic)
 Unpack the imported column indices and values, and combine into matrix. More...
 
std::string DistributorSendTypeEnumToString (EDistributorSendType sendType)
 Convert an EDistributorSendType enum value to a string. More...
 
std::string DistributorHowInitializedEnumToString (EDistributorHowInitialized how)
 Convert an EDistributorHowInitialized enum value to a string. More...
 
template<class SC , class LO , class GO , class NT >
void lclDotRaw (typename ::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *const resultRaw, const ::Tpetra::MultiVector< SC, LO, GO, NT > &X, const ::Tpetra::MultiVector< SC, LO, GO, NT > &Y, const bool resultOnDevice)
 Compute the local dot product(s), columnwise, of X and Y. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool isLocallyFitted (const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map1, const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &map2)
 Is map1 locally fitted to map2? More...
 
bool congruent (const Teuchos::Comm< int > &comm1, const Teuchos::Comm< int > &comm2)
 Whether the two communicators are congruent. More...
 
template<class DualViewType >
Teuchos::ArrayView< typename DualViewType::t_dev::value_type > getArrayViewFromDualView (const DualViewType &x)
 Get a Teuchos::ArrayView which views the host Kokkos::View of the input 1-D Kokkos::DualView. More...
 
template<class T , class DT >
Kokkos::DualView< T *, DT > getDualViewCopyFromArrayView (const Teuchos::ArrayView< const T > &x_av, const char label[], const bool leaveOnHost)
 Get a 1-D Kokkos::DualView which is a deep copy of the input Teuchos::ArrayView (which views host memory). More...
 
int countPackTriplesCount (const ::Teuchos::Comm< int > &comm, int &size, std::ostream *errStrm=NULL)
 Compute the buffer size required by packTriples for packing the number of matrix entries ("triples"). More...
 
int packTriplesCount (const int numEnt, char outBuf[], const int outBufSize, int &outBufCurPos, const ::Teuchos::Comm< int > &comm, std::ostream *errStrm=NULL)
 Pack the count (number) of matrix triples. More...
 
int unpackTriplesCount (const char inBuf[], const int inBufSize, int &inBufCurPos, int &numEnt, const ::Teuchos::Comm< int > &comm, std::ostream *errStrm=NULL)
 Unpack just the count of triples from the given input buffer. More...
 
template<class ScalarType , class OrdinalType >
int countPackTriples (const int numEnt, const ::Teuchos::Comm< int > &comm, int &size, std::ostream *errStrm=NULL)
 Compute the buffer size required by packTriples for packing numEnt number of (i,j,A(i,j)) matrix entries ("triples"). More...
 
template<class ScalarType , class OrdinalType >
int packTriples (const OrdinalType [], const OrdinalType [], const ScalarType [], const int, char [], const int, int &, const ::Teuchos::Comm< int > &, std::ostream *errStrm=NULL)
 Pack matrix entries ("triples" (i, j, A(i,j))) into the given output buffer. More...
 
template<class ScalarType , class OrdinalType >
int unpackTriples (const char [], const int, int &, OrdinalType [], OrdinalType [], ScalarType [], const int, const ::Teuchos::Comm< int > &, std::ostream *errStrm=NULL)
 Unpack matrix entries ("triples" (i, j, A(i,j))) from the given input buffer. More...
 
template<class SC , class GO >
int readAndDealOutTriples (std::istream &inputStream, std::size_t &curLineNum, std::size_t &totalNumEntRead, std::function< int(const GO, const GO, const SC &)> processTriple, const std::size_t maxNumEntPerMsg, const ::Teuchos::Comm< int > &comm, const bool tolerant=false, std::ostream *errStrm=NULL, const bool debug=false)
 On Process 0 in the given communicator, read sparse matrix entries (in chunks of at most maxNumEntPerMsg entries at a time) from the input stream, and "deal them out" to all other processes in the communicator. More...
 

Detailed Description

Namespace for Tpetra implementation details.

Warning
Do NOT rely on the contents of this namespace.

Enumeration Type Documentation

◆ EStorageStatus

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

When a CrsGraph or CrsMatrix is not fill complete, its data live in one of three storage formats:

  1. "2-D storage": The graph stores column indices as "array of arrays," and the matrix stores values as "array of arrays." The graph must have k_numRowEntries_ allocated. This only ever exists if the graph was created with DynamicProfile. A matrix with 2-D storage must own its graph, and the graph must have 2-D storage.

  2. "Unpacked 1-D storage": The graph uses a row offsets array, and stores column indices in a single array. The matrix also stores values in a single array. "Unpacked" means that there may be extra space in each row: that is, the row offsets array only says how much space there is in each row. The graph must use k_numRowEntries_ to find out how many entries there actually are in the row. A matrix with unpacked 1-D storage must own its graph, and the graph must have unpacked 1-D storage.

  3. "Packed 1-D storage": The matrix may or may not own the graph. "Packed" means that there is no extra space in each row. Thus, the k_numRowEntries_ array is not necessary and may have been deallocated. If the matrix was created with a constant ("static") graph, this must be true.

With respect to the Kokkos refactor version of Tpetra, "2-D storage" should be considered a legacy option.

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 storage is "not optimized," we retain the 1-D unpacked or 2-D format, and thus retain this enum value.

Definition at line 197 of file Tpetra_CrsGraph_decl.hpp.

◆ EDistributorSendType

The type of MPI send that Distributor should use.

This is an implementation detail of Distributor. Please do not rely on these values in your code.

Definition at line 76 of file Tpetra_Distributor.hpp.

◆ EDistributorHowInitialized

Enum indicating how and whether a Distributor was initialized.

This is an implementation detail of Distributor. Please do not rely on these values in your code.

Definition at line 94 of file Tpetra_Distributor.hpp.

Function Documentation

◆ convertColumnIndicesFromGlobalToLocal()

template<class LO , class GO , class DT , class OffsetType , class NumEntType >
OffsetType Tpetra::Details::convertColumnIndicesFromGlobalToLocal ( const Kokkos::View< LO *, DT > &  lclColInds,
const Kokkos::View< const GO *, DT > &  gblColInds,
const Kokkos::View< const OffsetType *, DT > &  ptr,
const LocalMap< LO, GO, DT > &  lclColMap,
const Kokkos::View< const NumEntType *, DT > &  numRowEnt 
)

Convert a (StaticProfile) CrsGraph's global column indices into local column indices.

Parameters
lclColInds[out] On output: The graph's local column indices. This may alias gblColInds, if LO == GO.
gblColInds[in] On input: The graph's global column indices. This may alias lclColInds, if LO == GO.
ptr[in] The graph's row offsets.
lclColMap[in] "Local" (threaded-kernel-worthy) version of the column Map.
numRowEnt[in] Array with number of entries in each row.
Returns
the number of "bad" global column indices (that don't live in the column Map on the calling process).

Definition at line 151 of file Tpetra_CrsGraph_def.hpp.

◆ computeOffsetsFromCounts()

template<class OffsetsViewType , class CountsViewType , class SizeType = typename OffsetsViewType::size_type>
OffsetsViewType::non_const_value_type Tpetra::Details::computeOffsetsFromCounts ( const OffsetsViewType &  ptr,
const CountsViewType &  counts 
)

Compute offsets from counts.

Compute offsets from counts via prefix sum:

ptr[i+1] = {j=0}^{i} counts[j]

Thus, ptr[i+1] - ptr[i] = counts[i], so that ptr[i+1] = ptr[i] + counts[i]. If we stored counts[i] in ptr[i+1] on input, then the formula is ptr[i+1] += ptr[i].

Returns
Sum of all counts; last entry of ptr.
Template Parameters
OffsetsViewTypeType of the Kokkos::View specialization used to store the offsets; the output array of this function.
CountsViewTypeType of the Kokkos::View specialization used to store the counts; the input array of this function.
SizeTypeThe parallel loop index type; a built-in integer type. Defaults to the type of the input View's dimension. You may use a shorter type to improve performance.

The type of each entry of the ptr array must be able to store the sum of all the entries of counts. This functor makes no attempt to check for overflow in this sum.

Definition at line 278 of file Tpetra_Details_computeOffsets.hpp.

◆ computeOffsetsFromConstantCount()

template<class OffsetsViewType , class CountType , class SizeType = typename OffsetsViewType::size_type>
OffsetsViewType::non_const_value_type Tpetra::Details::computeOffsetsFromConstantCount ( const OffsetsViewType &  ptr,
const CountType &  count 
)

Compute offsets from a constant count.

Compute offsets from a constant count via prefix sum:

ptr[i+1] = {j=0}^{i} count

Thus, ptr[i+1] - ptr[i] = count, so that ptr[i+1] = ptr[i] + count.

Returns
Sum of all counts; last entry of ptr.
Template Parameters
OffsetsViewTypeType of the Kokkos::View specialization used to store the offsets; the output array of this function.
CountTypeType of the constant count; the input argument of this function.
SizeTypeThe parallel loop index type; a built-in integer type. Defaults to the type of the output View's dimension. You may use a shorter type to improve performance.

The type of each entry of the ptr array must be able to store ptr.dimension_0 () * count. This functor makes no attempt to check for overflow in this sum.

Definition at line 407 of file Tpetra_Details_computeOffsets.hpp.

◆ copyConvert()

template<class OutputViewType , class InputViewType >
void Tpetra::Details::copyConvert ( const OutputViewType &  dst,
const InputViewType &  src 
)

Copy values from the 1-D Kokkos::View src, to the 1-D Kokkos::View dst, of the same length. The entries of src and dst may have different types, but it must be possible to copy-construct each entry of dst with its corresponding entry of src.

Everything above is an implementation detail of this function, copyConvert.

Definition at line 289 of file Tpetra_Details_copyConvert.hpp.

◆ copyOffsets()

template<class OutputViewType , class InputViewType >
void Tpetra::Details::copyOffsets ( const OutputViewType &  dst,
const InputViewType &  src 
)

Copy row offsets (in a sparse graph or matrix) from src to dst. The offsets may have different types.

The implementation reserves the right to do bounds checking if the offsets in the two arrays have different types.

Everything above is an implementation detail of this function, copyOffsets. This function in turn is an implementation detail of FixedHashTable, in particular of the "copy constructor" that copies a FixedHashTable from one Kokkos device to another. copyOffsets copies the array of offsets (ptr_).

Definition at line 407 of file Tpetra_Details_copyOffsets.hpp.

◆ create_mirror_view_from_raw_host_array()

template<class ValueType , class OutputDeviceType >
Impl::CreateMirrorViewFromUnmanagedHostArray<ValueType, OutputDeviceType>::output_view_type Tpetra::Details::create_mirror_view_from_raw_host_array ( const OutputDeviceType &  ,
ValueType *  inPtr,
const size_t  inSize,
const bool  copy = true,
const char  label[] = "" 
)

Variant of Kokkos::create_mirror_view that takes a raw host 1-d array as input.

Given a pointer to a 1-D array in host memory, and the number of entries in the array, return a Kokkos::View that lives in OutputDeviceType, and that is a mirror view of the input array. By default, copy the host data to the output View, if necessary.

Definition at line 201 of file Tpetra_Details_createMirrorView.hpp.

◆ crsMatrixSumIntoValues_sortedSortedLinear()

template<class SparseMatrixType , class ValsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type Tpetra::Details::crsMatrixSumIntoValues_sortedSortedLinear ( const SparseMatrixType &  A,
const typename SparseMatrixType::ordinal_type  lclRow,
const typename SparseMatrixType::ordinal_type  lclColInds[],
const typename SparseMatrixType::ordinal_type  sortPerm[],
const ValsViewType &  vals,
const typename SparseMatrixType::ordinal_type  numEntInInput,
const bool  forceAtomic = false,
const bool  checkInputIndices = true 
)

A(lclRow, lclColsInds[sortPerm[j]]) += vals[sortPerm[j]], for all j in 0 .. eltDim-1.

In the row of the matrix A with the local row index lclRow, find entries with column indices lclColInds, and sum into those entries with vals. Assume that lclColInds[sortPerm] is sorted, and that the column indices in that row of the matrix are sorted as well. Use linear search to find the entries in that row of the matrix.

Template Parameters
SparseMatrixTypeSpecialization of KokkosSparse::CrsMatrix.
ValsViewTypeSpecialization of a 1-D Kokkos::View.
Parameters
A[in/out] Sparse matrix whose entries to modify.
lclRow[in] Local index of the row in the matrix A to modify. lclRow MUST be a valid local row index of A.
lclColInds[in] Local column indices to modify in that row.
sortPerm[in] Permutation that makes lclColInds sorted. That is, lclColInds[sortPerm] is sorted.
vals[in] Input 1-D Kokkos::View of the values to ruse. This is a Kokkos::View and not a raw 1-D array, because it may be strided, if the original element being used (see crsMatrixSumInElement) has a column-major layout.
numEntInInput[in] Number of entries in the input. This function will read the first numEntInInput entries of lclColInds, sortPerm, and vals.
forceAtomic[in] Whether to use atomic updates when modifying the entries of the matrix A. This MUST be a compile-time constant. It defaults to whether the matrix's Kokkos execution space is NOT Kokkos::Serial.
checkInputIndices[in] Whether to check whether the input indices are valid column indices before just using them. For forwards compatibility, this should always be a compile-time constant. Default is true, that is, always check.
Returns
If checkInputIndices is true, return the number of input indices that are valid column indices in that row of the matrix. If checkInputIndices is false, just return numEntInInput.

Definition at line 93 of file Tpetra_Details_crsMatrixAssembleElement.hpp.

◆ crsMatrixReplaceValues_sortedSortedLinear()

template<class SparseMatrixType , class ValsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type Tpetra::Details::crsMatrixReplaceValues_sortedSortedLinear ( const SparseMatrixType &  A,
const typename SparseMatrixType::ordinal_type  lclRow,
const typename SparseMatrixType::ordinal_type  lclColInds[],
const typename SparseMatrixType::ordinal_type  sortPerm[],
const ValsViewType &  vals,
const typename SparseMatrixType::ordinal_type  numEntInInput,
const bool  forceAtomic = false,
const bool  checkInputIndices = true 
)

A(lclRow, lclColsInds[sortPerm[j]]) = vals[sortPerm[j]], for all j in 0 .. eltDim-1.

In the row of the matrix A with the local row index lclRow, find entries with column indices lclColInds, and replace those entries with vals. Assume that lclColInds[sortPerm] is sorted, and that the column indices in that row of the matrix are sorted as well. Use linear search to find the entries in that row of the matrix.

Template Parameters
SparseMatrixTypeSpecialization of KokkosSparse::CrsMatrix.
ValsViewTypeSpecialization of a 1-D Kokkos::View.
Parameters
A[in/out] Sparse matrix whose entries to modify.
lclRow[in] Local index of the row in the matrix A to modify. lclRow MUST be a valid local row index of A.
lclColInds[in] Local column indices to modify in that row.
sortPerm[in] Permutation that makes lclColInds sorted. That is, lclColInds[sortPerm] is sorted.
vals[in] Input 1-D Kokkos::View of the values to use. This is a Kokkos::View and not a raw 1-D array, because it may be strided, if the original element being used (see crsMatrixSumInElement) has a column-major layout.
numEntInInput[in] Number of entries in the input. This function will read the first numEntInInput entries of lclColInds, sortPerm, and vals.
forceAtomic[in] Whether to use atomic updates when modifying the entries of the matrix A. For forwards compatibility, this should always be a compile-time constant. It defaults to whether the matrix's Kokkos execution space is NOT Kokkos::Serial.
checkInputIndices[in] Whether to check whether the input indices are valid column indices before just using them. This MUST be a compile-time constant. Default is true, that is, always check.
Returns
If checkInputIndices is true, return the number of input indices that are valid column indices in that row of the matrix. If checkInputIndices is false, just return numEntInInput.

Definition at line 217 of file Tpetra_Details_crsMatrixAssembleElement.hpp.

◆ crsMatrixAssembleElement_sortedLinear()

template<class SparseMatrixType , class VectorViewType , class RhsViewType , class LhsViewType >
KOKKOS_FUNCTION SparseMatrixType::ordinal_type Tpetra::Details::crsMatrixAssembleElement_sortedLinear ( const SparseMatrixType &  A,
const VectorViewType &  x,
typename SparseMatrixType::ordinal_type  lids[],
typename SparseMatrixType::ordinal_type  sortPerm[],
const RhsViewType &  rhs,
const LhsViewType &  lhs,
const bool  forceAtomic = false,
const bool  checkInputIndices = true 
)

A(lids[j], lids[j]) += lhs(j,j) and x(lids[j]) += rhs(j), for all j in 0 .. eltDim-1.

Assume the following:

  • In each row of the sparse matrix A, the column indices are sorted.
  • The row and column indices of A have the same local indexing scheme, that is, a valid row index is a valid column index and vice versa.

Sum the dense "element" matrix (2-D Kokkos::View) lhs into the entries of the sparse matrix A corresponding to the input row and column indices lids. Also, sum the dense "element" vector (1-D Kokkos::View) rhs into the entries of the dense vector x corresponding to the input row indices lids.

Template Parameters
SparseMatrixTypeSpecialization of KokkosSparse::CrsMatrix.
RhsViewTypeSpecialization of a 1-D Kokkos::View.
LhsViewTypeSpecialization of a 2-D Kokkos::View.
Parameters
A[in/out] Sparse matrix (KokkosSparse::CrsMatrix) to modify.
x[in/out] Dense vector (1-D Kokkos::View) to modify.
lids[in/out] Local row and column indices of A to modify. This function may sort this array, and output the permutation that makes it sorted to sortPerm. lids must have the same number of entries as rhs.dimension_0(), lhs.dimension_0(), and lhs.dimension_1().
sortPerm[out] Permutation that makes lids (on input) sorted. It must have the same number of writeable entries as lids (see above).
rhs[in] Dense "element" vector of input values to sum into the dense vector x; a 1-D Kokkos::View. It must have the same number of entries as each dimension of lhs.
lhs[in] Dense, square "element" matrix of input values to sum into the sparse matrix A; a 2-D Kokkos::View. Each of its dimensions must be the same as the number of entries in rhs.
forceAtomic[in] Whether to use atomic updates when modifying the entries of the matrix A and vector x. For forwards compatibility, this should always be a compile-time constant. It defaults to whether the matrix's Kokkos execution space is NOT Kokkos::Serial.
checkInputIndices[in] Whether to check whether the input indices are valid column indices before just using them. This MUST be a compile-time constant. Default is true, that is, always check.
Returns
If checkInputIndices is true, return the number of input indices that are valid column indices in that row of the matrix. If checkInputIndices is false, just return numEntInInput.

Definition at line 356 of file Tpetra_Details_crsMatrixAssembleElement.hpp.

◆ gathervPrint()

void Tpetra::Details::gathervPrint ( std::ostream &  out,
const std::string &  s,
const Teuchos::Comm< int > &  comm 
)

On Process 0 in the given communicator, print strings from each process in that communicator, in rank order.

For each process in the given communicator comm, send its string s to Process 0 in that communicator. Process 0 prints the strings in rank order.

This is a collective over the given communicator comm. Process 0 promises not to store all the strings in its memory. This function's total memory usage on any process is proportional to the calling process' string length, plus the max string length over any process. This does NOT depend on the number of processes in the communicator. Thus, we call this a "memory-scalable" operation. While the function's name suggests MPI_Gatherv, the implementation may NOT use MPI_Gather or MPI_Gatherv, because neither of those are not memory scalable.

Process 0 prints nothing other than what is in the string. It does not add an endline after each string, nor does it identify each string with its owning process' rank. If you want either of those in the string, you have to put it there yourself.

Parameters
out[out] The output stream to which to write. ONLY Process 0 in the given communicator will write to this. Thus, this stream need only be valid on Process 0.
s[in] The string to write. Each process in the given communicator has its own string. Strings may be different on different processes. Zero-length strings are OK.
comm[in] The communicator over which this operation is a collective.

Definition at line 52 of file Tpetra_Details_gathervPrint.cpp.

◆ getDiagCopyWithoutOffsets()

template<class DiagType , class LocalMapType , class CrsMatrixType >
static LocalMapType::local_ordinal_type Tpetra::Details::getDiagCopyWithoutOffsets ( const DiagType &  D,
const LocalMapType &  rowMap,
const LocalMapType &  colMap,
const CrsMatrixType &  A 
)
static

Given a locally indexed, local sparse matrix, and corresponding local row and column Maps, extract the matrix's diagonal entries into a 1-D Kokkos::View.

Warning
This is an implementation detail of Tpetra::CrsMatrix. This function may disappear or change its interface at any time.

This function implements much of the one-argument overload of Tpetra::CrsMatrix::getLocalDiagCopy, for the case where the matrix is fill complete. The function computes offsets of diagonal entries inline, and does not store them. If you want to store the offsets, call computeOffsets() instead.

Template Parameters
DiagType1-D nonconst Kokkos::View
CrsMatrixTypeSpecialization of KokkosSparse::CrsMatrix
LocalMapTypeSpecialization of Tpetra::Details::LocalMap; type of the "local" part of a Tpetra::Map
Parameters
D[out] 1-D Kokkos::View to which to write the diagonal entries.
rowMap[in] "Local" part of the sparse matrix's row Map.
colMap[in] "Local" part of the sparse matrix's column Map.
A[in] The sparse matrix.

Definition at line 182 of file Tpetra_Details_getDiagCopyWithoutOffsets_decl.hpp.

◆ getLocalDiagCopyWithoutOffsetsNotFillComplete()

template<class SC , class LO , class GO , class NT >
LO Tpetra::Details::getLocalDiagCopyWithoutOffsetsNotFillComplete ( ::Tpetra::Vector< SC, LO, GO, NT > &  diag,
const ::Tpetra::RowMatrix< SC, LO, GO, NT > &  A,
const bool  debug = false 
)

Given a locally indexed, global sparse matrix, extract the matrix's diagonal entries into a Tpetra::Vector.

Warning
This is an implementation detail of Tpetra::CrsMatrix. This function may disappear or change its interface at any time.

This function is a work-around for Github Issue #499. It implements one-argument Tpetra::CrsMatrix::getLocalDiagCopy for the case where the matrix is not fill complete. The function computes offsets of diagonal entries inline, and does not store them. If you want to store the offsets, call computeOffsets() instead.

Template Parameters
SCSame as first template parameter (Scalar) of Tpetra::CrsMatrix and Tpetra::Vector.
LOSame as second template parameter (LocalOrdinal) of Tpetra::CrsMatrix and Tpetra::Vector.
GOSame as third template parameter (GlobalOrdinal) of Tpetra::CrsMatrix and Tpetra::Vector.
NTSame as fourth template parameter (Node) of Tpetra::CrsMatrix and Tpetra::Vector.
Parameters
diag[out] Tpetra::Vector to which to write the diagonal entries. Its Map must be the same (in the sense of Tpetra::Map::isSameAs()) as the row Map of A.
A[in] The sparse matrix. Must be a Tpetra::RowMatrix (the base class of Tpetra::CrsMatrix), must be locally indexed, and must have row views.
debug[in] Whether to do extra run-time checks. This costs MPI communication. The default is false in a release build, and true in a debug build.

We pass in the sparse matrix as a Tpetra::RowMatrix because the implementation of Tpetra::CrsMatrix uses this function, and we want to avoid a circular header dependency. On the other hand, the implementation does not actually depend on Tpetra::CrsMatrix.

Definition at line 192 of file Tpetra_Details_getDiagCopyWithoutOffsets_def.hpp.

◆ iallreduce()

template<class InputViewType , class OutputViewType >
std::shared_ptr<CommRequest> Tpetra::Details::iallreduce ( const InputViewType &  sendbuf,
const OutputViewType &  recvbuf,
const ::Teuchos::EReductionType  op,
const ::Teuchos::Comm< int > &  comm 
)

Nonblocking all-reduce, for either rank-1 or rank-0 Kokkos::View objects.

Template Parameters
InputViewTypeType of the send buffer
OutputViewTypeType of the receive buffer

This function wraps MPI_Iallreduce. It does a nonblocking all-reduce over the input communicator comm, from sendbuf into recvbuf, using op as the reduction operator. The function returns without blocking; the all-reduce only blocks for completion when one calls wait() on the returned request.

Parameters
sendbuf[in] Input buffer; must be either a rank-1 or rank-0 Kokkos::View, and must have the same rank as recvbuf.
recvbuf[in] Output buffer; must be either a rank-1 or rank-0 Kokkos::View, and must have the same rank as sendbuf.
op[in] Teuchos enum representing the reduction operator.
comm[in] Communicator over which to do the all-reduce.

sendbuf and recvbuf must either be disjoint, or identical (point to the same array). They may not partially overlap. Furthermore, if they are identical, the input communicator must be an intracommunicator. It may not be an intercommunicator. If you don't know what an intercommunicator is, you probably just have an intracommunicator, so everything is fine.

Definition at line 813 of file Tpetra_Details_iallreduce.hpp.

◆ isInterComm()

bool Tpetra::Details::isInterComm ( const Teuchos::Comm< int > &  comm)

Return true if and only if the input communicator wraps an MPI intercommunicator.

The most common MPI communicators are intracommunicators ("in<i>tra</i>," not "in<i>ter</i>"). This includes MPI_COMM_WORLD, MPI_COMM_SELF, and the results of MPI_Comm_dup and MPI_Comm_split. Intercommunicators come from MPI_Intercomm_create.

This distinction matters because only collectives over intracommunicators may use MPI_IN_PLACE, to let the send and receive buffers alias each other. Collectives over intercommunicators may not use MPI_IN_PLACE.

Parameters
comm[in] The input communicator.
Returns
Whether the input communicator wraps an MPI intercommunicator.

Definition at line 66 of file Tpetra_Details_isInterComm.cpp.

◆ makeColMap()

template<class LO , class GO , class NT >
int Tpetra::Details::makeColMap ( Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &  colMap,
const Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &  domMap,
const RowGraph< LO, GO, NT > &  graph,
const bool  sortEachProcsGids = true,
std::ostream *  errStrm = NULL 
)

Make the graph's column Map.

Template Parameters
LOLocal ordinal type; the type of local indices in the graph.
GOGlobal ordinal type; the type of global indices in the graph.
NTNode type; the third template parameter of Tpetra::CrsGraph.
Parameters
colMap[out] On output: pointer to the column Map for the given graph. This is only valid if the returned error code is zero on all processes in the communicator of domMap (see below). This may be the same (literally the same object as) as domMap, depending on the graph.
domMap[in] The domain Map to use for creating the column Map. This need not be the same as graph.getDomainMap(). It's OK for the latter to be null, in fact. domMap needs to be passed in by RCP, because it's possible for the returned column Map colMap (see above) to equal domMap.
graph[in] The graph for which to make a column Map. This function does NOT modify the graph's column Map, if it happens to have one already. Thus, this function supports graph modification.
sortEachProcsGids[in] Whether to sort column Map GIDs associated with each remote process in ascending order. This is true by default. If false, leave GIDs in their original order as discovered in the graph by iterating in ascending order through the local rows of the graph.
errStrm[out] If nonnull, print error messages to this.
Returns
Error code; zero if and only if successful. This value is local to the calling process. On error, the code may be zero on some processes, and nonzero on other processes. You, the user, are responsible for propagating that error state to all processes.

This function always makes a column Map, even if the Map already has one. This makes it possible to change the graph's structure, and have its column Map and corresponding Import update in the same way.

The sortEachProcsGids argument corresponds to sortGhostsAssociatedWithEachProcessor_ in CrsGraph. This function always groups remote GIDs by process rank, so that all remote GIDs with the same owning rank occur contiguously. The sortEachProcsGids argument (see above) whether this function sorts remote GIDs in increasing order within those groups. This function sorts by default. This behavior differs from Epetra, which does not sort remote GIDs with the same owning process. means "sort remote GIDs." If you don't want to sort, for compatibility with Epetra, set sortEachProcsGids to false.

Definition at line 66 of file Tpetra_Details_makeColMap_def.hpp.

◆ makeOptimizedColMap()

template<class MapType >
MapType Tpetra::Details::makeOptimizedColMap ( std::ostream &  errStream,
bool &  lclErr,
const MapType &  domMap,
const MapType &  colMap 
)

Return an optimized reordering of the given column Map.

Template Parameters
MapTypeA specialization of Map.
Parameters
err[out] Output stream for human-readable error reporting. This is local to the calling process and may differ on different processes.
lclErr[out] On output: true if anything went wrong on the calling process. This value is local to the calling process and may differ on different processes.
domMap[in] Domain Map of a CrsGraph or CrsMatrix.
colMap[in] Original column Map of the same CrsGraph or CrsMatrix as domMap.
Returns
The possibly reordered column Map newColMap.

This is a convenience wrapper for makeOptimizedColMapAndImport(). (Please refer to that function's documentation in this file.) It does everything that that function does, except that it does not compute a new Import.

Definition at line 310 of file Tpetra_Details_makeOptimizedColMap.hpp.

◆ makeOptimizedColMapAndImport()

template<class MapType >
std::pair<MapType, Teuchos::RCP<typename OptColMap<MapType>::import_type> > Tpetra::Details::makeOptimizedColMapAndImport ( std::ostream &  errStream,
bool &  lclErr,
const MapType &  domMap,
const MapType &  colMap,
const typename OptColMap< MapType >::import_type *  oldImport,
const bool  makeImport 
)

Return an optimized reordering of the given column Map. Optionally, recompute an Import from the input domain Map to the new column Map.

Template Parameters
MapTypeA specialization of Map.

This function takes a domain Map and a column Map of a distributed graph (Tpetra::CrsGraph) or matrix (e.g., Tpetra::CrsMatrix). It then creates a new column Map, which optimizes the performance of an Import operation from the domain Map to the new column Map. This function also optionally creates that Import. Creating the new column Map and its Import at the same time saves some communication, since making the Import requires some of the same information that optimizing the column Map does.

Parameters
err[out] Output stream for human-readable error reporting. This is local to the calling process and may differ on different processes.
lclErr[out] On output: true if anything went wrong on the calling process. This value is local to the calling process and may differ on different processes.
domMap[in] Domain Map of a CrsGraph or CrsMatrix.
colMap[in] Original column Map of the same CrsGraph or CrsMatrix as domMap.
makeImport[in] Whether to make and return an Import from the input domain Map to the new column Map.
Returns
The possibly reordered column Map newColMap, and the corresponding Import from domMap to newColMap. The latter is nonnull if and only if makeImport is true.
Precondition
domMap and colMap must have the same or congruent communicators.
On all calling processes, the indices in colMap must be a subset of the indices in domMap.

The returned column Map's global indices (GIDs) will have the following order on all calling processes:

  1. GIDs that occur in both colMap and domMap (on the calling process) go first.
  2. GIDs in colMap on the calling process, but not in the domain Map on the calling process, follow. They are ordered first contiguously by their owning process rank (in the domain Map), then in increasing order within that.

This imitates the ordering used by AztecOO and Epetra. Storing indices contiguously that are owned by the same process (in the domain Map) permits the use of contiguous send and receive buffers in Distributor, which is used in an Import operation.

Definition at line 383 of file Tpetra_Details_makeOptimizedColMap.hpp.

◆ countMergeUnsortedIndices()

template<class OrdinalType , class IndexType >
IndexType Tpetra::Details::countMergeUnsortedIndices ( const OrdinalType  curInds[],
const IndexType  numCurInds,
const OrdinalType  inputInds[],
const IndexType  numInputInds 
)

Count the number of column indices that can be merged into the current row, assuming that both the current row's indices and the input indices are unsorted.

Neither the current row's entries, nor the input, are sorted. Return the number of input entries that can be merged into the current row. Don't actually merge them. 'numCurInds' corresponds to 'midPos' in mergeUnsortedIndices.

The current indices are NOT allowed to have repeats, but the input indices ARE allowed to have repeats. (The whole point of these methods is to keep the current entries without repeats – "merged in.") Repeats in the input are counted separately with respect to merges.

The unsorted case is bad for asymptotics, but the asymptotics only show up with dense or nearly dense rows, which are bad for other reasons.

Definition at line 74 of file Tpetra_Details_Merge.hpp.

◆ countMergeSortedIndices()

template<class OrdinalType , class IndexType >
IndexType Tpetra::Details::countMergeSortedIndices ( const OrdinalType  curInds[],
const IndexType  numCurInds,
const OrdinalType  inputInds[],
const IndexType  numInputInds 
)

Count the number of column indices that can be merged into the current row, assuming that both the current row's indices and the input indices are sorted.

Both the current row's entries and the input are sorted. Return the number of input entries that can be merged into the current row. Don't actually merge them. 'numCurInds' corresponds to 'midPos' in mergeSortedIndices.

The current indices are NOT allowed to have repeats, but the input indices ARE allowed to have repeats. (The whole point of these methods is to keep the current entries without repeats – "merged in.") Repeats in the input are counted separately with respect to merges.

The sorted case is good for asymptotics, but imposes an order on the entries of each row. Sometimes users don't want that.

Definition at line 133 of file Tpetra_Details_Merge.hpp.

◆ mergeSortedIndices()

template<class OrdinalType , class IndexType >
std::pair<bool, IndexType> Tpetra::Details::mergeSortedIndices ( OrdinalType  curInds[],
const IndexType  midPos,
const IndexType  endPos,
const OrdinalType  inputInds[],
const IndexType  numInputInds 
)

Attempt to merge the input indices into the current row's column indices, assuming that both the current row's indices and the input indices are sorted.

Both the current row's entries and the input are sorted. If and only if the current row has enough space for the input (after merging), merge the input with the current row.

Assume that both curInds and inputInds are sorted. Current indices: curInds[0 .. midPos-1]. Extra space at end: curInds[midPos .. endPos-1] Input indices to merge in: inputInds[0 .. numInputInds]. Any of those could be empty.

If the merge succeeded, return true and the new number of entries in the row. Else, return false and the new number of entries in the row required to fit the input.

The sorted case is good for asymptotics, but imposes an order on the entries of each row. Sometimes users don't want that.

Definition at line 205 of file Tpetra_Details_Merge.hpp.

◆ mergeUnsortedIndices()

template<class OrdinalType , class IndexType >
std::pair<bool, IndexType> Tpetra::Details::mergeUnsortedIndices ( OrdinalType  curInds[],
const IndexType  midPos,
const IndexType  endPos,
const OrdinalType  inputInds[],
const IndexType  numInputInds 
)

Attempt to merge the input indices into the current row's column indices, assuming that both the current row's indices and the input indices are unsorted.

Neither the current row's entries nor the input are sorted. If and only if the current row has enough space for the input (after merging), merge the input with the current row.

Assume that neither curInds nor inputInds are sorted. Current indices: curInds[0 .. midPos-1]. Extra space at end: curInds[midPos .. endPos-1] Input indices to merge in: inputInds[0 .. numInputInds]. Any of those could be empty.

If the merge succeeded, return true and the new number of entries in the row. Else, return false and the new number of entries in the row required to fit the input.

The unsorted case is bad for asymptotics, but the asymptotics only show up with dense or nearly dense rows, which are bad for other reasons.

Definition at line 325 of file Tpetra_Details_Merge.hpp.

◆ mergeUnsortedIndicesAndValues()

template<class OrdinalType , class ValueType , class IndexType >
std::pair<bool, IndexType> Tpetra::Details::mergeUnsortedIndicesAndValues ( OrdinalType  curInds[],
ValueType  curVals[],
const IndexType  midPos,
const IndexType  endPos,
const OrdinalType  inputInds[],
const ValueType  inputVals[],
const IndexType  numInputInds 
)

Attempt to merge the input indices and values into the current row's column indices and corresponding values, assuming that both the current row's indices and the input indices are unsorted.

Neither the current row's entries nor the input are sorted. If and only if the current row has enough space for the input (after merging), merge the input with the current row.

Assume that neither curInds nor inputInds are sorted. Current indices: curInds[0 .. midPos-1]. Current values: curVals[0 .. midPos-1]. Extra space for indices at end: curInds[midPos .. endPos-1]. Extra space for values at end: curVals[midPos .. endPos-1]. Input indices to merge in: inputInds[0 .. numInputInds]. Input values to merge in: inputVals[0 .. numInputInds].

If the merge succeeded, return true and the new number of entries in the row. Else, return false and the new number of entries in the row required to fit the input.

The unsorted case is bad for asymptotics, but the asymptotics only show up with dense or nearly dense rows, which are bad for other reasons.

Definition at line 419 of file Tpetra_Details_Merge.hpp.

◆ packCrsMatrixRow()

template<class LocalMatrixType , class LocalMapType >
KOKKOS_FUNCTION bool Tpetra::Details::packCrsMatrixRow ( const LocalMatrixType &  lclMatrix,
const LocalMapType &  lclColMap,
char *const  numEntOut,
char *const  valOut,
char *const  indOut,
const size_t  numEnt,
const typename LocalMatrixType::ordinal_type  lclRow 
)

Packs a single row of the CrsMatrix.

Data (bytes) describing the row of the CRS matrix are "packed" (concatenated) in to a single char* as

LO number of entries | GO column indices > – number of entries | column indices | values – SC values |

Template Parameters
LocalMatrixTypethe specialization of the KokkosSparse::CrsMatrix local matrix
LocalMapTypethe type of the local column map

Definition at line 526 of file Tpetra_Details_packCrsMatrix.hpp.

◆ packCrsMatrix()

template<class LocalMatrixType , class LocalMapType >
bool Tpetra::Details::packCrsMatrix ( const LocalMatrixType &  lclMatrix,
const LocalMapType &  lclColMap,
std::unique_ptr< std::string > &  errStr,
Teuchos::Array< char > &  exports,
const Teuchos::ArrayView< size_t > &  numPacketsPerLID,
size_t &  constantNumPackets,
const Teuchos::ArrayView< const typename LocalMatrixType::ordinal_type > &  exportLIDs,
const int  myRank,
Distributor  
)

Pack specified entries of the given local sparse matrix for communication.

Warning
This is an implementation detail of Tpetra::CrsMatrix.
Parameters
errStr[in/out] If an error occurs on any participating process, allocate this if it is null, then fill the string with local error reporting. This is purely local to the process that reports the error. The caller is responsible for synchronizing across processes.
exports[in/out] Output pack buffer; resized if needed.
numPacketsPerLID[out] Entry k gives the number of bytes packed for row exportLIDs[k] of the local matrix.
exportLIDs[in] Local indices of the rows to pack.
lclMatrix[in] The local sparse matrix to pack.
Returns
true if no errors occurred on the calling process, else false. This is purely local to the process that discovered the error. The caller is responsible for synchronizing across processes.

Definition at line 632 of file Tpetra_Details_packCrsMatrix.hpp.

◆ shortSortKeysAndValues_2()

template<class KeyType , class ValueType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeysAndValues_2 ( KeyType  keys[2],
ValueType  values[2] 
)

Sort keys and values jointly, by keys, for arrays of length 2.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
ValueTypeCopy constructible, assignable
Parameters
keys[in/out] Length 2 array of keys. This function sorts this keys array, and applies the same permutation to the values array.
values[in/out] Length 2 array of values.

Definition at line 138 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeys_2()

template<class KeyType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeys_2 ( KeyType  keys[2])

Sort length-2 array of keys.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
Parameters
keys[in/out] Length-2 array of keys to sort.

Definition at line 154 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeysAndValues_3()

template<class KeyType , class ValueType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeysAndValues_3 ( KeyType  keys[3],
ValueType  values[3] 
)

Sort keys and values jointly, by keys, for arrays of length 3.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
ValueTypeCopy constructible, assignable
Parameters
keys[in/out] Length 3 array of keys. This function sorts this keys array, and applies the same permutation to the values array.
values[in/out] Length 3 array of values.

Definition at line 174 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeys_3()

template<class KeyType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeys_3 ( KeyType  keys[3])

Sort length-3 array of keys.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
Parameters
keys[in/out] Length-3 array of keys to sort.

Definition at line 196 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeysAndValues_4()

template<class KeyType , class ValueType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeysAndValues_4 ( KeyType  keys[4],
ValueType  values[4] 
)

Sort keys and values jointly, by keys, for arrays of length 4.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
ValueTypeCopy constructible, assignable
Parameters
keys[in/out] Length 4 array of keys. This function sorts this keys array, and applies the same permutation to the values array.
values[in/out] Length 4 array of values.

Definition at line 222 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeys_4()

template<class KeyType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeys_4 ( KeyType  keys[4])

Sort length-4 array of keys.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
Parameters
keys[in/out] Length-4 array of keys to sort.

Definition at line 246 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeysAndValues_8()

template<class KeyType , class ValueType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeysAndValues_8 ( KeyType  keys[8],
ValueType  values[8] 
)

Sort keys and values jointly, by keys, for arrays of length 8.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
ValueTypeCopy constructible, assignable
Parameters
keys[in/out] Length 8 array of keys. This function sorts this keys array, and applies the same permutation to the values array.
values[in/out] Length 8 array of values.

Definition at line 274 of file Tpetra_Details_shortSort.hpp.

◆ shortSortKeys_8()

template<class KeyType >
KOKKOS_FUNCTION void Tpetra::Details::shortSortKeys_8 ( KeyType  keys[8])

Sort length-8 array of keys.

Template Parameters
KeyTypeGreater-than comparable, copy constructible, assignable
Parameters
keys[in/out] Length-8 array of keys to sort.

Definition at line 312 of file Tpetra_Details_shortSort.hpp.

◆ shellSortKeysAndValues()

template<class KeyType , class ValueType , class IndexType >
KOKKOS_FUNCTION void Tpetra::Details::shellSortKeysAndValues ( KeyType  keys[],
ValueType  values[],
const IndexType  n 
)

Shellsort (yes, it's one word) the input array keys, and apply the resulting permutation to the input array values.

mfh 28 Nov 2016, 17 Dec 2016: I adapted this function from sh_sort2 in Tpetra_Util.hpp (in this directory).

Definition at line 349 of file Tpetra_Details_shortSort.hpp.

◆ shellSortKeys()

template<class KeyType , class IndexType >
KOKKOS_FUNCTION void Tpetra::Details::shellSortKeys ( KeyType  keys[],
const IndexType  n 
)

Shellsort (yes, it's one word) the input array keys.

Parameters
keys[in/out] Input array of keys to sort.
n[in] Length of the input array keys.

Definition at line 390 of file Tpetra_Details_shortSort.hpp.

◆ unpackCrsMatrixAndCombine()

template<class LocalMatrixType , class LocalMapType >
bool Tpetra::Details::unpackCrsMatrixAndCombine ( LocalMatrixType &  lclMatrix,
const LocalMapType &  lclColMap,
std::unique_ptr< std::string > &  errStr,
const Teuchos::ArrayView< const typename LocalMatrixType::ordinal_type > &  importLIDs,
const Teuchos::ArrayView< const char > &  imports,
const Teuchos::ArrayView< const size_t > &  numPacketsPerLID,
size_t  constantNumPackets,
const int  myRank,
Distributor ,
CombineMode  combineMode,
const bool  atomic 
)

Unpack the imported column indices and values, and combine into matrix.

Warning
The allowed combineMode are: ADD, REPLACE, and ABSMAX. INSERT is not allowed.

Definition at line 336 of file Tpetra_Details_unpackCrsMatrix.hpp.

◆ DistributorSendTypeEnumToString()

std::string Tpetra::Details::DistributorSendTypeEnumToString ( EDistributorSendType  sendType)

Convert an EDistributorSendType enum value to a string.

This is an implementation detail of Distributor. Please do not rely on this function in your code.

Definition at line 49 of file Tpetra_Distributor.cpp.

◆ DistributorHowInitializedEnumToString()

std::string Tpetra::Details::DistributorHowInitializedEnumToString ( EDistributorHowInitialized  how)

Convert an EDistributorHowInitialized enum value to a string.

This is an implementation detail of Distributor. Please do not rely on this function in your code.

Definition at line 70 of file Tpetra_Distributor.cpp.

◆ lclDotRaw()

template<class SC , class LO , class GO , class NT >
void Tpetra::Details::lclDotRaw ( typename ::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *const  resultRaw,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  X,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  Y,
const bool  resultOnDevice 
)

Compute the local dot product(s), columnwise, of X and Y.

Warning
This is an implementation detail of Tpetra::idot. Users should never call this function.

This implements the following cases:

  • If X and Y have the same number of columns (zero or more), compute the dot products of corresponding columns of X and Y. That is, resultRaw[j] = dot(X(:,j), Y(:,j)).
  • If X has one column and Y has more than one column, compute the dot products of X and each column of Y in turn. That is, resultRaw[j] = dot(X(:,0), Y(:,j)).
  • If X has more than one column and Y has one column, compute the dot products of each column of X in turn with X. That is, resultRaw[j] = dot(X(:,j), Y(:,0)).
Template Parameters
SCSame as the first template parameter of Tpetra::MultiVector.
LOSame as the second template parameter of Tpetra::MultiVector.
GOSame as the third template parameter of Tpetra::MultiVector.
NTSame as the fourth template parameter of Tpetra::MultiVector.
Parameters
resultRaw[out] Raw pointer to output array of dot products.
X[in] First input MultiVector.
Y[in] Second input MultiVector.
resultOnDevice[in] Whether resultRaw points to memory accessible from device. If not, it may only be accessed from a host execution space.

Definition at line 103 of file Tpetra_idot.hpp.

◆ isLocallyFitted()

template<class LocalOrdinal , class GlobalOrdinal , class Node >
bool Tpetra::Details::isLocallyFitted ( const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &  map1,
const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &  map2 
)

Is map1 locally fitted to map2?

Parameters
map1[in] The first Map
map2[in] The second Map

For Map instances map1 and map2, we say that map1 is locally fitted to map2 (on the calling process), when the initial indices of map1 (on the calling process) are the same and in the same order as those of map2. "Fittedness" is entirely a local (per MPI process) property.

The predicate "is map1 fitted to map2 ?" is not symmetric. For example, map2 may have more entries than map1.

Fittedness on a process can let Tpetra avoid deep copies in some Export or Import (communication) operations. Tpetra could use this, for example, in optimizing its sparse matrix-vector multiply.

Definition at line 1925 of file Tpetra_Map_def.hpp.

◆ congruent()

bool Tpetra::Details::congruent ( const Teuchos::Comm< int > &  comm1,
const Teuchos::Comm< int > &  comm2 
)

Whether the two communicators are congruent.

Two communicators are congruent when they have the same number of processes, and those processes occur in the same rank order.

If both communicators are MpiComm instances, this function returns true exactly when MPI_Comm_compare returns MPI_IDENT (the communicators are handles for the same object) or MPI_CONGRUENT. SerialComm instances are always congruent. An MpiComm is congruent to a SerialComm if the MpiComm has only one process. This function is symmetric in its arguments.

If either Comm instance is neither an MpiComm nor a SerialComm, this method cannot do any better than to compare their process counts.

Two communicators are congruent when they have the same number of processes, and those processes occur in the same rank order.

If both communicators are Teuchos::MpiComm instances, this function returns true exactly when MPI_Comm_compare returns MPI_IDENT (the communicators are handles for the same object) or MPI_CONGRUENT on their MPI_Comm handles. Any two Teuchos::SerialComm instances are always congruent. An MpiComm instance is congruent to a SerialComm instance if and only if the MpiComm has one process. This function is symmetric in its arguments.

If either Teuchos::Comm instance is neither an MpiComm nor a SerialComm, this method cannot do any better than to compare their process counts.

Definition at line 65 of file Tpetra_Util.cpp.

◆ getArrayViewFromDualView()

template<class DualViewType >
Teuchos::ArrayView<typename DualViewType::t_dev::value_type> Tpetra::Details::getArrayViewFromDualView ( const DualViewType &  x)

Get a Teuchos::ArrayView which views the host Kokkos::View of the input 1-D Kokkos::DualView.

Precondition
The input DualView must be sync'd to host.
Parameters
x[in] A specialization of Kokkos::DualView.
Returns
Teuchos::ArrayView that views the host version of the DualView's data.

Definition at line 875 of file Tpetra_Util.hpp.

◆ getDualViewCopyFromArrayView()

template<class T , class DT >
Kokkos::DualView<T*, DT> Tpetra::Details::getDualViewCopyFromArrayView ( const Teuchos::ArrayView< const T > &  x_av,
const char  label[],
const bool  leaveOnHost 
)

Get a 1-D Kokkos::DualView which is a deep copy of the input Teuchos::ArrayView (which views host memory).

Template Parameters
TThe type of the entries of the input Teuchos::ArrayView.
DTThe Kokkos Device type.
Parameters
x_av[in] The Teuchos::ArrayView to copy.
label[in] String label for the Kokkos::DualView.
leaveOnHost[in] If true, the host version of the returned Kokkos::DualView is most recently updated (and the DualView may need a sync to device). If false, the device version is most recently updated (and the DualView may need a sync to host).
Returns
Kokkos::DualView that is a deep copy of the input Teuchos::ArrayView.

Definition at line 909 of file Tpetra_Util.hpp.

◆ countPackTriplesCount()

int Tpetra::Details::countPackTriplesCount ( const ::Teuchos::Comm< int > &  comm,
int &  size,
std::ostream *  errStrm = NULL 
)

Compute the buffer size required by packTriples for packing the number of matrix entries ("triples").

countPackTriples tells me an upper bound on how much buffer space I need to hold numEnt triples. packTriplesCount actually packs numEnt, the number of triples. countPackTriplesCount tells me an upper bound on how much buffer space I need to hold the number of triples, not the triples themselves.

Parameters
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
size[out] Pack buffer size in bytes (sizeof(char)).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 173 of file Tpetra_Details_PackTriples.cpp.

◆ packTriplesCount()

int Tpetra::Details::packTriplesCount ( const int  numEnt,
char  outBuf[],
const int  outBufSize,
int &  outBufCurPos,
const ::Teuchos::Comm< int > &  comm,
std::ostream *  errStrm = NULL 
)

Pack the count (number) of matrix triples.

This function is NOT the same thing as countPackTriples. countPackTriples tells me an upper bound on how much buffer space I need to hold numEnt triples. packTriplesCount actually packs numEnt, the number of triples. countPackTriplesCount tells me an upper bound on how much buffer space I need to hold the number of triples, not the triples themselves.

Parameters
numEnt[in] Number of matrix entries ("triples") to pack.
outBuf[out] Output buffer.
outBufSize[out] Total output buffer size in bytes.
outBufCurPos[in/out] Current position from which to start writing to the output buffer. This corresponds to the 'position' in/out argument of MPI_Pack.
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 201 of file Tpetra_Details_PackTriples.cpp.

◆ unpackTriplesCount()

int Tpetra::Details::unpackTriplesCount ( const char  inBuf[],
const int  inBufSize,
int &  inBufCurPos,
int &  numEnt,
const ::Teuchos::Comm< int > &  comm,
std::ostream *  errStrm = NULL 
)

Unpack just the count of triples from the given input buffer.

We store the count of triples as an int, because MPI buffer sizes are int.

Parameters
inBuf[in] Input buffer.
inBufSize[out] Total input buffer size in bytes.
inBufCurPos[in/out] Current position from which to start reading from the input buffer. This corresponds to the 'position' in/out argument of MPI_Unpack.
numEnt[out] Number of matrix entries ("triples") that were packed.
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 232 of file Tpetra_Details_PackTriples.cpp.

◆ countPackTriples()

template<class ScalarType , class OrdinalType >
int Tpetra::Details::countPackTriples ( const int  numEnt,
const ::Teuchos::Comm< int > &  comm,
int &  size,
std::ostream *  errStrm = NULL 
)

Compute the buffer size required by packTriples for packing numEnt number of (i,j,A(i,j)) matrix entries ("triples").

This function is NOT the same thing as packTriplesCount. countPackTriples tells me an upper bound on how much buffer space I need to hold numEnt triples. packTriplesCount actually packs numEnt, the number of triples. countPackTriplesCount tells me an upper bound on how much buffer space I need to hold the number of triples, not the triples themselves.

Template Parameters
ScalarTypeType of each matrix entry A(i,j).
OrdinalTypeType of each matrix index i or j.
Parameters
numEnt[in] Number of matrix entries ("triples") to pack.
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
size[out] Pack buffer size in bytes (sizeof(char)).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 353 of file Tpetra_Details_PackTriples.hpp.

◆ packTriples()

template<class ScalarType , class OrdinalType >
int Tpetra::Details::packTriples ( const OrdinalType  [],
const OrdinalType  [],
const ScalarType  [],
const int  ,
char  [],
const int  ,
int &  ,
const ::Teuchos::Comm< int > &  ,
std::ostream *  errStrm = NULL 
)

Pack matrix entries ("triples" (i, j, A(i,j))) into the given output buffer.

Template Parameters
ScalarTypeType of each matrix entry A(i,j).
OrdinalTypeType of each matrix index i or j.
Parameters
gblRowInds[in] Row indices to pack.
gblColInds[in] Column indices to pack.
val[in] Matrix values A(i,j) to pack.
numEnt[in] Number of matrix entries ("triples") to pack.
outBuf[out] Output buffer.
outBufSize[out] Total output buffer size in bytes.
outBufCurPos[in/out] Current position from which to start writing to the output buffer. This corresponds to the 'position' in/out argument of MPI_Pack.
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 463 of file Tpetra_Details_PackTriples.hpp.

◆ unpackTriples()

template<class ScalarType , class OrdinalType >
int Tpetra::Details::unpackTriples ( const char  [],
const int  ,
int &  ,
OrdinalType  [],
OrdinalType  [],
ScalarType  [],
const int  ,
const ::Teuchos::Comm< int > &  ,
std::ostream *  errStrm = NULL 
)

Unpack matrix entries ("triples" (i, j, A(i,j))) from the given input buffer.

Template Parameters
ScalarTypeType of each matrix entry A(i,j).
OrdinalTypeType of each matrix index i or j.
Parameters
inBuf[in] Input buffer.
inBufSize[out] Total pack buffer size in bytes (sizeof(char)).
inBufCurPos[in/out] Current position from which to start reading from the input buffer. This corresponds to the 'position' in/out argument of MPI_Unpack.
gblRowInds[out] Row indices unpacked.
gblColInds[out] Column indices unpacked.
val[out] Matrix values A(i,j) unpacked.
numEnt[in] Number of matrix entries ("triples") to unpack. If you don't know it, then you should have senders pack the triples count as the first thing in the buffer, and unpack it first via unpackTriplesCount().
comm[in] Communicator used in sending and receiving the packed entries. (MPI wants this, so we have to include it.).
errStrm[out] If nonnull, print any error messages to this stream, else don't print error messages.
Returns
Error code. MPI_SUCCESS (0) if successful, else nonzero.
Warning
It only makes sense to call this function if using MPI. If not building with MPI, this function is a stub that returns nonzero.

Definition at line 579 of file Tpetra_Details_PackTriples.hpp.

◆ readAndDealOutTriples()

template<class SC , class GO >
int Tpetra::Details::readAndDealOutTriples ( std::istream &  inputStream,
std::size_t &  curLineNum,
std::size_t &  totalNumEntRead,
std::function< int(const GO, const GO, const SC &)>  processTriple,
const std::size_t  maxNumEntPerMsg,
const ::Teuchos::Comm< int > &  comm,
const bool  tolerant = false,
std::ostream *  errStrm = NULL,
const bool  debug = false 
)

On Process 0 in the given communicator, read sparse matrix entries (in chunks of at most maxNumEntPerMsg entries at a time) from the input stream, and "deal them out" to all other processes in the communicator.

This is a collective over the communicator.

Template Parameters
SCThe type of the value of each matrix entry.
GOThe type of each (global) index of each matrix entry.
Parameters
inputStream[in/out] Input stream from which to read Matrix Market - format matrix entries ("triples"). Only Process 0 in the communicator needs to be able to access this.
curLineNum[in/out] On both input and output, the current line number in the input stream. (In the Matrix Market format, sparse matrix entries cannot start until at least line 3 of the file.) This is only valid on Process 0.
totalNumEntRead[out] Total number of matrix entries (triples) read on Process 0. This is only valid on Process 0.
processTriple[in] Closure, generally with side effects, that takes in and stores off a sparse matrix entry. First argument is the (global) row index, second argument is the (global) column index, and third argument is the value of the entry. The closure must NOT do MPI communication. Return value is an error code, that is zero if and only if the closure succeeded. We intend for you to use this to call CooMatrix::insertEntry.
comm[in] Communicator to use for receiving the triples.
tolerant[in] Whether to read tolerantly.
errStrm[in] If not NULL, print any error messages to this stream.
Returns
Error code; 0 if and only if success.

Definition at line 939 of file Tpetra_Details_ReadTriples.hpp.