Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Namespaces | Classes | Typedefs | Enumerations | Functions
Tpetra Namespace Reference

Namespace Tpetra contains the class and methods constituting the Tpetra library. More...

Namespaces

 Details
 Nonmember function that computes a residual Computes R = B - A * X.
 
 Impl
 Namespace for new Tpetra features that are not ready for public release, but are ready for evaluation by friendly expert users.
 
 Ext
 Namespace for external Tpetra functionality.
 
 MatrixMatrix
 Distributed sparse matrix-matrix multiply and add.
 
 TripleMatrixMultiply
 Distributed sparse triple matrix product.
 

Classes

class  BlockCrsMatrix
 Sparse matrix whose entries are small dense square blocks, all of the same dimensions. More...
 
class  BlockMultiVector
 MultiVector for multiple degrees of freedom per mesh point. More...
 
class  BlockVector
 Vector for multiple degrees of freedom per mesh point. More...
 
class  project1st
 Binary function that returns its first argument. More...
 
class  project2nd
 Binary function that returns its second argument. More...
 
class  ScopeGuard
 Scope guard whose destructor automatically calls Tpetra::finalize for you. More...
 
struct  RowInfo
 Allocation information for a locally owned row in a CrsGraph or CrsMatrix. More...
 
class  CrsGraph
 A distributed graph accessed by rows (adjacency lists) and stored sparsely. More...
 
class  CrsMatrix
 Sparse matrix that presents a row-oriented interface that lets users read or modify entries. More...
 
class  CrsMatrixMultiplyOp
 A class for wrapping a CrsMatrix multiply in a Operator. More...
 
class  Directory
 Implement mapping from global ID to process ID and local ID. More...
 
class  DistObject
 Base class for distributed Tpetra objects that support data redistribution. More...
 
class  Distributor
 Sets up and executes a communication plan for a Tpetra DistObject. More...
 
class  Export
 Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distribution. More...
 
class  FECrsGraph
 A distributed graph accessed by rows (adjacency lists) and stored sparsely. More...
 
class  Import
 Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distribution. More...
 
class  ImportExportData
 Implementation detail of Import and Export. More...
 
class  LocalCrsMatrixOperator
 Abstract interface for local operators (e.g., matrices and preconditioners). More...
 
class  LocalOperator
 Abstract interface for local operators (e.g., matrices and preconditioners). More...
 
class  Map
 A parallel distribution of indices over processes. More...
 
class  MultiVector
 One or more distributed dense vectors. More...
 
class  Operator
 Abstract interface for operators (e.g., matrices and preconditioners). More...
 
class  Packable
 Abstract base class for objects that can be the source of an Import or Export operation, and that also know how to pack their data to send to the target object. More...
 
class  RowGraph
 An abstract interface for graphs accessed by rows. More...
 
class  RowMatrix
 A read-only, row-oriented interface to a sparse matrix. More...
 
class  RowMatrixTransposer
 Construct and (optionally) redistribute the explicitly stored transpose of a CrsMatrix. More...
 
class  SrcDistObject
 Abstract base class for objects that can be the source of an Import or Export operation. More...
 
class  Vector
 A distributed dense vector. More...
 
class  CrsMatrixStruct
 Struct that holds views of the contents of a CrsMatrix. More...
 

Typedefs

typedef Teuchos_Ordinal Array_size_type
 Size type for Teuchos Array objects. More...
 
typedef size_t global_size_t
 Global size_t object. More...
 
template<class LocalAccessType >
using with_local_access_function_argument_type = typename Details::GetNonowningLocalObject< LocalAccessType >::nonowning_local_object_type
 Type of the local object, that is an argument to the function the user gives to withLocalAccess. More...
 

Enumerations

enum  CombineMode {
  ADD, INSERT, REPLACE, ABSMAX,
  ZERO
}
 Rule for combining data in an Import or Export. More...
 
enum  LocalGlobal
 Enum for local versus global allocation of Map entries. More...
 
enum  LookupStatus { AllIDsPresent, IDNotPresent }
 Return status of Map remote index lookup (getRemoteIndexList()). More...
 
enum  ProfileType
 
enum  OptimizeOption { DoOptimizeStorage, DoNotOptimizeStorage }
 
enum  ESweepDirection
 Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR). More...
 
enum  EScaling
 Whether "scaling" a matrix means multiplying or dividing its entries. More...
 

Functions

template<class CrsMatrixType >
void applyDirichletBoundaryConditionToLocalMatrixRows (const typename CrsMatrixType::execution_space &execSpace, CrsMatrixType &A, const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, typename CrsMatrixType::device_type > &lclRowInds)
 For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on device. More...
 
template<class CrsMatrixType >
void applyDirichletBoundaryConditionToLocalMatrixRows (CrsMatrixType &A, const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, typename CrsMatrixType::device_type > &lclRowInds)
 For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on device, using the default execution space instance. More...
 
template<class CrsMatrixType >
void applyDirichletBoundaryConditionToLocalMatrixRows (CrsMatrixType &A, const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, Kokkos::HostSpace > &lclRowInds)
 For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on host, using the default host execution space instance. More...
 
template<class Scalar , class LO , class GO , class Node >
void blockCrsMatrixWriter (BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::string const &fileName)
 Helper function to write a BlockCrsMatrix. Calls the 3-argument version. More...
 
template<class Scalar , class LO , class GO , class Node >
void blockCrsMatrixWriter (BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::string const &fileName, Teuchos::ParameterList const &params)
 Helper function to write a BlockCrsMatrix. Calls the 3-argument version. More...
 
template<class Scalar , class LO , class GO , class Node >
void blockCrsMatrixWriter (BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os)
 Helper function to write a BlockCrsMatrix. Calls the 3-argument version. More...
 
template<class Scalar , class LO , class GO , class Node >
void blockCrsMatrixWriter (BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os, Teuchos::ParameterList const &params)
 Helper function to write a BlockCrsMatrix. More...
 
template<class Scalar , class LO , class GO , class Node >
void writeMatrixStrip (BlockCrsMatrix< Scalar, LO, GO, Node > const &A, std::ostream &os, Teuchos::ParameterList const &params)
 Helper function called by blockCrsMatrixWriter. More...
 
template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< BlockCrsMatrix
< Scalar, LO, GO, Node > > 
convertToBlockCrsMatrix (const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &pointMatrix, const LO &blockSize)
 Non-member constructor that creates a BlockCrsMatrix from an existing point CrsMatrix. More...
 
template<class LO , class GO , class Node >
Teuchos::RCP< const
Tpetra::Map< LO, GO, Node > > 
createMeshMap (LO const &blockSize, const Tpetra::Map< LO, GO, Node > &pointMap)
 Helper function to generate a mesh map from a point map. Important! It's assumed that point GIDs associated with a single mesh GID appear consecutively in pointMap. More...
 
template<class ViewType , class CoefficientType , class LayoutType = typename ViewType::array_layout, class IndexType = int, const int rank = ViewType::rank>
KOKKOS_INLINE_FUNCTION void SCAL (const CoefficientType &alpha, const ViewType &x)
 x := alpha*x, where x is either rank 1 (a vector) or rank 2 (a matrix). More...
 
template<class ViewType , class InputType , class LayoutType = typename ViewType::array_layout, class IndexType = int, const int rank = ViewType::rank>
KOKKOS_INLINE_FUNCTION void FILL (const ViewType &x, const InputType &val)
 Set every entry of x to val. More...
 
template<class CoefficientType , class ViewType1 , class ViewType2 , class LayoutType1 = typename ViewType1::array_layout, class LayoutType2 = typename ViewType2::array_layout, class IndexType = int, const int rank = ViewType1::rank>
KOKKOS_INLINE_FUNCTION void AXPY (const CoefficientType &alpha, const ViewType1 &x, const ViewType2 &y)
 y := y + alpha * x (dense vector or matrix update) More...
 
template<class ViewType1 , class ViewType2 , class LayoutType1 = typename ViewType1::array_layout, class LayoutType2 = typename ViewType2::array_layout, class IndexType = int, const int rank = ViewType1::rank>
KOKKOS_INLINE_FUNCTION void COPY (const ViewType1 &x, const ViewType2 &y)
 Deep copy x into y, where x and y are either rank 1 (vectors) or rank 2 (matrices) with the same dimension(s). More...
 
template<class VecType1 , class BlkType , class VecType2 , class CoeffType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void GEMV (const CoeffType &alpha, const BlkType &A, const VecType1 &x, const VecType2 &y)
 y := y + alpha * A * x (dense matrix-vector multiply) More...
 
template<class ViewType1 , class ViewType2 , class ViewType3 , class CoefficientType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void GEMM (const char transA[], const char transB[], const CoefficientType &alpha, const ViewType1 &A, const ViewType2 &B, const CoefficientType &beta, const ViewType3 &C)
 Small dense matrix-matrix multiply: C := alpha*A*B + beta*C More...
 
template<class LittleBlockType , class LittleVectorType >
KOKKOS_INLINE_FUNCTION void GETF2 (const LittleBlockType &A, const LittleVectorType &ipiv, int &info)
 Computes A = P*L*U. More...
 
template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void GETRS (const char mode[], const LittleBlockType &A, const LittleIntVectorType &ipiv, const LittleScalarVectorType &B, int &info)
 Solve the linear system(s) AX=B, using the result of GETRF or GETF2. More...
 
template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void GETRI (const LittleBlockType &A, const LittleIntVectorType &ipiv, const LittleScalarVectorType &work, int &info)
 Compute inverse of A, using result of GETRF or GETF2. More...
 
void setCombineModeParameter (Teuchos::ParameterList &plist, const std::string &paramName)
 Set CombineMode parameter in a Teuchos::ParameterList. More...
 
std::string combineModeToString (const CombineMode combineMode)
 Human-readable string representation of the given CombineMode. More...
 
template<class SC , class LO , class GO , class NT >
Details::EquilibrationInfo
< typename Kokkos::ArithTraits
< SC >::val_type, typename
NT::device_type > 
computeRowOneNorms (const Tpetra::RowMatrix< SC, LO, GO, NT > &A)
 Compute global row one-norms ("row sums") of the input sparse matrix A, in a way suitable for one-sided (left only) equilibration. More...
 
template<class SC , class LO , class GO , class NT >
Details::EquilibrationInfo
< typename Kokkos::ArithTraits
< SC >::val_type, typename
NT::device_type > 
computeRowAndColumnOneNorms (const Tpetra::RowMatrix< SC, LO, GO, NT > &A, const bool assumeSymmetric)
 Compute global row and column one-norms ("row sums" and "column sums") of the input sparse matrix A, in a way suitable for two-sided (left and right) equilibration. More...
 
bool isInitialized ()
 Whether Tpetra is in an initialized state. More...
 
Teuchos::RCP< const
Teuchos::Comm< int > > 
getDefaultComm ()
 Get Tpetra's default communicator. More...
 
void initialize (int *argc, char ***argv)
 Initialize Tpetra. More...
 
void initialize (int *argc, char ***argv, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Initialize Tpetra. More...
 
void finalize ()
 Finalize Tpetra. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsGraph
< LocalOrdinal, GlobalOrdinal,
Node > > 
createCrsGraph (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &map, size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember function to create an empty CrsGraph given a row Map and the max number of entries allowed locally per row. More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &importer, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsGraph constructor that fuses Import and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > importAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &rowImporter, const Import< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &domainImporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsGraph constructor that fuses Import and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &exporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsGraph constructor that fuses Export and fillComplete(). More...
 
template<class CrsGraphType >
Teuchos::RCP< CrsGraphType > exportAndFillCompleteCrsGraph (const Teuchos::RCP< const CrsGraphType > &sourceGraph, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &rowExporter, const Export< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > &domainExporter, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsGraphType::local_ordinal_type, typename CrsGraphType::global_ordinal_type, typename CrsGraphType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsGraph constructor that fuses Export and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &importer, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsMatrix constructor that fuses Import and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > importAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &rowImporter, const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainImporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsMatrix constructor that fuses Import and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &exporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap=Teuchos::null, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Nonmember CrsMatrix constructor that fuses Export and fillComplete(). More...
 
template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > exportAndFillCompleteCrsMatrix (const Teuchos::RCP< const CrsMatrixType > &sourceMatrix, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &rowExporter, const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &domainExporter, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &domainMap, const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > &params)
 Nonmember CrsMatrix constructor that fuses Export and fillComplete(). More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
createCrsMatrix (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t maxNumEntriesPerRow=0, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
 Create an empty CrsMatrix given a row map and a single integer upper bound on the number of stored entries per row. More...
 
template<class OpScalar , class MatScalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP
< CrsMatrixMultiplyOp
< OpScalar, MatScalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
createCrsMatrixMultiplyOp (const Teuchos::RCP< const CrsMatrix< MatScalar, LocalOrdinal, GlobalOrdinal, Node > > &A)
 Non-member function to create a CrsMatrixMultiplyOp. More...
 
template<class DistObjectType >
void removeEmptyProcessesInPlace (Teuchos::RCP< DistObjectType > &input, const Teuchos::RCP< const Map< typename DistObjectType::local_ordinal_type, typename DistObjectType::global_ordinal_type, typename DistObjectType::node_type > > &newMap)
 Remove processes which contain no elements in this object's Map. More...
 
template<class DistObjectType >
void removeEmptyProcessesInPlace (Teuchos::RCP< DistObjectType > &input)
 Remove processes which contain no elements in this object's Map. More...
 
Teuchos::Array< std::string > distributorSendTypes ()
 Valid values for Distributor's "Send type" parameter. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Export
< LocalOrdinal, GlobalOrdinal,
Node > > 
createExport (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &src, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &tgt)
 Nonmember "constructor" for Export objects. More...
 
template<class ExecutionSpace , class GlobalDataStructure , class UserFunctionType >
void for_each (const char kernelLabel[], ExecutionSpace execSpace, GlobalDataStructure &X, UserFunctionType f)
 Apply a function entrywise to each local entry of a Tpetra global data structure, analogously to std::for_each. More...
 
template<class GlobalDataStructure , class UserFunctionType >
void for_each (const char kernelLabel[], GlobalDataStructure &X, UserFunctionType f)
 Overload of for_each (see above) that runs on X's default Kokkos execution space. More...
 
template<class SC , class LO , class GO , class NT >
std::shared_ptr
< ::Tpetra::Details::CommRequest
idot (typename::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *resultRaw, const ::Tpetra::MultiVector< SC, LO, GO, NT > &X, const ::Tpetra::MultiVector< SC, LO, GO, NT > &Y)
 Nonblocking dot product, with either Tpetra::MultiVector or Tpetra::Vector inputs, and raw pointer or raw array output. More...
 
template<class SC , class LO , class GO , class NT >
std::shared_ptr
< ::Tpetra::Details::CommRequest
idot (const Kokkos::View< typename::Tpetra::Vector< SC, LO, GO, NT >::dot_type, typename::Tpetra::Vector< SC, LO, GO, NT >::device_type > &result, const ::Tpetra::Vector< SC, LO, GO, NT > &X, const ::Tpetra::Vector< SC, LO, GO, NT > &Y)
 Nonblocking dot product, with Tpetra::Vector inputs, and rank-0 (single value) Kokkos::View output. More...
 
template<class SC , class LO , class GO , class NT >
std::shared_ptr
< ::Tpetra::Details::CommRequest
idot (const Kokkos::View< typename::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *, typename::Tpetra::MultiVector< SC, LO, GO, NT >::device_type > &result, const ::Tpetra::MultiVector< SC, LO, GO, NT > &X, const ::Tpetra::MultiVector< SC, LO, GO, NT > &Y)
 Nonblocking dot product, with Tpetra::MultiVector inputs, and rank-1 (one-dimensional array) Kokkos::View output. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
createImport (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &src, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &tgt)
 Nonmember constructor for Import. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Import
< LocalOrdinal, GlobalOrdinal,
Node > > 
createImport (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &src, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &tgt, const Teuchos::RCP< Teuchos::ParameterList > &plist)
 Nonmember constructor for Import that takes a ParameterList. More...
 
template<class SC , class LO , class GO , class NT >
void leftAndOrRightScaleCrsMatrix (Tpetra::CrsMatrix< SC, LO, GO, NT > &A, const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &rowScalingFactors, const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &colScalingFactors, const bool leftScale, const bool rightScale, const bool assumeSymmetric, const EScaling scaling)
 Left-scale and/or right-scale (in that order) the entries of the input Tpetra::CrsMatrix A. More...
 
template<class SC , class LO , class GO , class NT >
void leftAndOrRightScaleCrsMatrix (Tpetra::CrsMatrix< SC, LO, GO, NT > &A, const Tpetra::Vector< typename Kokkos::ArithTraits< SC >::mag_type, LO, GO, NT > &rowScalingFactors, const Tpetra::Vector< typename Kokkos::ArithTraits< SC >::mag_type, LO, GO, NT > &colScalingFactors, const bool leftScale, const bool rightScale, const bool assumeSymmetric, const EScaling scaling)
 Left-scale and/or right-scale (in that order) the entries of the input Tpetra::CrsMatrix A. More...
 
template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal > > 
createLocalMap (const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Nonmember constructor for a locally replicated Map with the default Kokkos Node. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createLocalMapWithNode (const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Nonmember constructor for a locally replicated Map with a specified Kokkos Node. More...
 
template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal > > 
createUniformContigMap (const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Non-member constructor for a uniformly distributed, contiguous Map with the default Kokkos Node. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createUniformContigMapWithNode (const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Non-member constructor for a uniformly distributed, contiguous Map with a user-specified Kokkos Node. More...
 
template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal > > 
createContigMap (const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Non-member constructor for a (potentially) non-uniformly distributed, contiguous Map using the default Kokkos::Device. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createContigMapWithNode (const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Nonmember constructor for a (potentially) nonuniformly distributed, contiguous Map for a user-specified, possibly nondefault Kokkos Node type. More...
 
template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal > > 
createNonContigMap (const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Nonmember constructor for a non-contiguous Map using the default Kokkos::Device type. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createNonContigMapWithNode (const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
 Nonmember constructor for a noncontiguous Map with a user-specified, possibly nondefault Kokkos Node type. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createOneToOne (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &M)
 Nonmember constructor for a contiguous Map with user-defined weights and a user-specified, possibly nondefault Kokkos Node type. More...
 
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map
< LocalOrdinal, GlobalOrdinal,
Node > > 
createOneToOne (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &M, const ::Tpetra::Details::TieBreak< LocalOrdinal, GlobalOrdinal > &tie_break)
 Creates a one-to-one version of the given Map where each GID lives on only one process. The given TieBreak object specifies the rule to break ties. More...
 
template<class DS , class DL , class DG , class DN , class SS , class SL , class SG , class SN >
void deep_copy (MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
 Copy the contents of the MultiVector src into dst. More...
 
template<class ST , class LO , class GO , class NT >
MultiVector< ST, LO, GO, NT > createCopy (const MultiVector< ST, LO, GO, NT > &src)
 Return a deep copy of the given MultiVector. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
createMultiVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, const size_t numVectors)
 Nonmember MultiVector "constructor": Create a MultiVector from a given Map. More...
 
template<class ST , class LO , class GO , class NT >
void deep_copy (MultiVector< ST, LO, GO, NT > &dst, const MultiVector< ST, LO, GO, NT > &src)
 Specialization of deep_copy for MultiVector objects with the same template parameters. More...
 
template<class SC , class LO , class GO , class NT >
LO replaceDiagonalCrsMatrix (::Tpetra::CrsMatrix< SC, LO, GO, NT > &matrix, const ::Tpetra::Vector< SC, LO, GO, NT > &newDiag)
 Replace diagonal entries of an input Tpetra::CrsMatrix matrix with values given in newDiag. More...
 
template<class ExecutionSpace , class GlobalDataStructure , class UnaryFunctionType >
void transform (const char kernelLabel[], ExecutionSpace execSpace, GlobalDataStructure &input, GlobalDataStructure &output, UnaryFunctionType f)
 For each local entry input_i of input, assign f(input_i) to the corresponding local entry output_i of output, analogously to std::transform. More...
 
template<class GlobalDataStructure , class UnaryFunctionType >
void transform (const char kernelLabel[], GlobalDataStructure &input, GlobalDataStructure &output, UnaryFunctionType f)
 Overload of transform (see above) that runs on the output object's default Kokkos execution space. More...
 
template<class ExecutionSpace , class GlobalDataStructure , class BinaryFunctionType >
void transform (const char kernelLabel[], ExecutionSpace execSpace, GlobalDataStructure &input1, GlobalDataStructure &input2, GlobalDataStructure &output, BinaryFunctionType f)
 Binary transform on the given execution space execSpace: output_i = f(input1_i, input2_i). More...
 
template<typename MapType , typename KeyArgType , typename ValueArgType >
MapType::iterator efficientAddOrUpdate (MapType &m, const KeyArgType &k, const ValueArgType &v)
 Efficiently insert or replace an entry in an std::map. More...
 
template<class IT1 , class IT2 >
void sort2 (const IT1 &first1, const IT1 &last1, const IT2 &first2)
 Sort the first array, and apply the resulting permutation to the second array. More...
 
template<class IT1 , class IT2 , class IT3 >
void sort3 (const IT1 &first1, const IT1 &last1, const IT2 &first2, const IT3 &first3)
 Sort the first array, and apply the same permutation to the second and third arrays. More...
 
template<class IT1 , class IT2 >
void merge2 (IT1 &indResultOut, IT2 &valResultOut, IT1 indBeg, IT1 indEnd, IT2 valBeg, IT2)
 Merge values in place, additively, with the same index. More...
 
template<class IT1 , class IT2 , class BinaryFunction >
void merge2 (IT1 &indResultOut, IT2 &valResultOut, IT1 indBeg, IT1 indEnd, IT2 valBeg, IT2 valEnd, BinaryFunction f)
 Merge values in place with the same index, using any associative binary function. More...
 
template<class KeyInputIterType , class ValueInputIterType , class KeyOutputIterType , class ValueOutputIterType , class BinaryFunction >
void keyValueMerge (KeyInputIterType keyBeg1, KeyInputIterType keyEnd1, ValueInputIterType valBeg1, ValueInputIterType valEnd1, KeyInputIterType keyBeg2, KeyInputIterType keyEnd2, ValueInputIterType valBeg2, ValueInputIterType valEnd2, KeyOutputIterType keyOut, ValueOutputIterType valOut, BinaryFunction f)
 Merge two sorted (by keys) sequences of unique (key,value) pairs by combining pairs with equal keys. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Vector< Scalar, LocalOrdinal,
GlobalOrdinal, Node > 
createCopy (const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &src)
 Return a deep copy of the given Vector. More...
 
template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Vector< Scalar,
LocalOrdinal, GlobalOrdinal,
Node > > 
createVector (const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map)
 Nonmember Vector "constructor": Create a Vector from a given Map. More...
 
template<class GlobalObjectType >
Details::LocalAccess
< GlobalObjectType,
Details::read_only
readOnly (GlobalObjectType &)
 Declare that you want to access the given global object's local data in read-only mode. More...
 
template<class GlobalObjectType >
Details::LocalAccess
< GlobalObjectType,
Details::read_only
readOnly (const GlobalObjectType &)
 Declare that you want to access the given global object's local data in read-only mode (overload for const GlobalObjectType). More...
 
template<class GlobalObjectType >
Details::LocalAccess
< GlobalObjectType,
Details::write_only
writeOnly (GlobalObjectType &)
 Declare that you want to access the given global object's local data in write-only mode. More...
 
template<class GlobalObjectType >
Details::LocalAccess
< GlobalObjectType,
Details::read_write
readWrite (GlobalObjectType &)
 Declare that you want to access the given global object's local data in read-and-write mode. More...
 
template<class... LocalAccessTypes>
void withLocalAccess (typename Details::ArgsToFunction< LocalAccessTypes...>::type userFunction, LocalAccessTypes...localAccesses)
 Get access to a Tpetra global object's local data. More...
 
Methods performing multiple matrix-vector products at once
template<class MatrixArray , class MultiVectorArray >
void batchedApply (const MatrixArray &Matrices, const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &X, MultiVectorArray &Y, typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type alpha=Teuchos::ScalarTraits< typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type >::one(), typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type beta=Teuchos::ScalarTraits< typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type >::zero(), Teuchos::RCP< Teuchos::ParameterList > params=Teuchos::null)
 Does multiply matrix apply() calls with a single X vector. More...
 

Detailed Description

Namespace Tpetra contains the class and methods constituting the Tpetra library.

Namespace for Tpetra classes and methods.

Typedef Documentation

typedef Teuchos_Ordinal Tpetra::Array_size_type

Size type for Teuchos Array objects.

Definition at line 51 of file Tpetra_ConfigDefs.hpp.

typedef size_t Tpetra::global_size_t

Global size_t object.

This type is intended to support scenarios where the global memory allocation is larger than that of a single node.

Currently, it is typedefed to size_t.

Definition at line 107 of file Tpetra_ConfigDefs.hpp.

template<class LocalAccessType >
using Tpetra::with_local_access_function_argument_type = typedef typename Details::GetNonowningLocalObject<LocalAccessType>:: nonowning_local_object_type

Type of the local object, that is an argument to the function the user gives to withLocalAccess.

Template Parameters
LocalAccessTypeSpecialization of LocalAccess.

withLocalAccess takes zero or more global objects, each wrapped in a LocalAccess struct. Wrapping happens indirectly, through e.g., the readOnly, writeOnly, or readWrite functions. withLocalAccess also takes a function, functor, or lambda with the same number of arguments (the "local objects") as the number of global objects. If users don't have C++14 (generic lambdas), then they will need to know the types of the local objects. This alias maps directly from the LocalAccess struct type, to the corresponding local object type.

Definition at line 248 of file Tpetra_withLocalAccess.hpp.

Enumeration Type Documentation

Rule for combining data in an Import or Export.

Import or Export (data redistribution) operations might need to combine data received from other processes with existing data on the calling process. This enum tells Tpetra how to do that for a specific Import or Export operation. Each Tpetra object may interpret the CombineMode in a different way, so you should check the Tpetra object's documentation for details.

Here is the list of supported combine modes:

  • ADD: Sum new values into existing values
  • INSERT: Insert new values that don't currently exist
  • REPLACE: Replace existing values with new values
  • ABSMAX: If $x_{old}$ is the old value and $x_{new}$ the incoming new value, replace $x_{old}$ with $\max\{ x_{old}, x_{new} \}$.
  • ZERO: Replace old values with zero

ADD and REPLACE are intended for modifying values that already exist. Tpetra objects will generally work correctly if those values don't already exist. (For example, ADD will behave like INSERT if the entry does not yet exist on the calling process.) However, performance may suffer.

The ZERO combine mode is a special case that bypasses communication. It may seem odd to include a "combine mode" that doesn't actually combine. However, this is useful for computations like domain decomposition with overlap. A ZERO combine mode with overlap is different than an ADD combine mode without overlap. (See Ifpack2::AdditiveSchwarz, which inspired inclusion of this combine mode.) Furthermore, Import and Export also encapsulate a local permutation; if you want only to execute the local permutation without communication, you may use the ZERO combine mode.

Enumerator
ADD 

Sum new values into existing values.

INSERT 

Insert new values that don't currently exist.

REPLACE 

Replace existing values with new values.

ABSMAX 

Replace old value with maximum of magnitudes of old and new values.

ZERO 

Replace old values with zero.

Definition at line 94 of file Tpetra_CombineMode.hpp.

Enum for local versus global allocation of Map entries.

LocallyReplicated means that the Map's entries are locally replicated across all processes.

GloballyDistributed means that the Map's entries are globally distributed across all processes.

Definition at line 116 of file Tpetra_ConfigDefs.hpp.

Return status of Map remote index lookup (getRemoteIndexList()).

Enumerator
AllIDsPresent 

All queried indices were present in the Map

IDNotPresent 

At least one of the specified indices was not present in the Map

Definition at line 122 of file Tpetra_ConfigDefs.hpp.

Allocation profile for matrix/graph entries

Definition at line 128 of file Tpetra_ConfigDefs.hpp.

Optimize storage option

Enumerator
DoOptimizeStorage 

Indicates that storage should be optimized

DoNotOptimizeStorage 

Indicates that storage should not be optimized

Definition at line 134 of file Tpetra_ConfigDefs.hpp.

Sweep direction for Gauss-Seidel or Successive Over-Relaxation (SOR).

Definition at line 237 of file Tpetra_ConfigDefs.hpp.

Whether "scaling" a matrix means multiplying or dividing its entries.

leftAndOrRightScaleCrsMatrix (see below) takes this as input.

Definition at line 60 of file Tpetra_leftAndOrRightScaleCrsMatrix_decl.hpp.

Function Documentation

template<class MatrixArray , class MultiVectorArray >
void Tpetra::batchedApply ( const MatrixArray &  Matrices,
const typename std::remove_pointer< typename MultiVectorArray::value_type >::type &  X,
MultiVectorArray &  Y,
typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type  alpha = Teuchos::ScalarTraits< typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::one(),
typename std::remove_pointer< typename MatrixArray::value_type >::type::scalar_type  beta = Teuchos::ScalarTraits< typename std::remove_pointer<typename MatrixArray::value_type>::type::scalar_type>::zero(),
Teuchos::RCP< Teuchos::ParameterList >  params = Teuchos::null 
)

Does multiply matrix apply() calls with a single X vector.

Computes Y[i] = beta*Y[i] + alpha * Matrices[i] * X, for i=1,...,N

This routine only communicates the interprocessor portions of X once, if such a thing is possible (aka the ColMap's of the matrices match).

If the Y's are replicated this will not use the reduced communication path.

If X is aliased to any Y but the last one, this routine will throw

Parameters
Matrices[in] - [std::vector|Teuchos::Array|Teuchos::ArrayRCP] of Tpetra::CrsMatrix objects. These matrices can different numbers of rows.
X[in] - Tpetra::MultiVector or Tpetra::Vector object.
Y[out] - [std::vector|Teuchos::Array|Teuchos::ArrayRCP] of Tpetra::MultiVector or Tpetra::Vector objects. These must have the same number of vectors as X.
alpha[in] - alpha parameter. Defaults to one.
beta[in] - beta parameter. Defaults to zero.
params[in/out] - The "can batch" parameter can either be unset, be true or be false on input. If it is unset, maps will be checked with isSameAs() for compatibility. If true, the map check will be skipped and batching will be used if none of the cheap checks fail. If false, batching will not be used. This parameter will be set on output to either true or false depending on if batching was used during this call. Defaults to NULL.

Definition at line 74 of file Tpetra_Apply_Helpers.hpp.

template<class CrsMatrixType >
void Tpetra::applyDirichletBoundaryConditionToLocalMatrixRows ( const typename CrsMatrixType::execution_space &  execSpace,
CrsMatrixType &  A,
const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, typename CrsMatrixType::device_type > &  lclRowInds 
)

For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on device.

Template Parameters
CrsMatrixTypeSpecialization of Tpetra::CrsMatrix.
Parameters
execSpace[in] Execution space instance on which to run.
A[in/out] Sparse matrix to modify.
lclRowInds[in] Local indices of the rows of A to modify.

Definition at line 197 of file Tpetra_applyDirichletBoundaryCondition.hpp.

template<class CrsMatrixType >
void Tpetra::applyDirichletBoundaryConditionToLocalMatrixRows ( CrsMatrixType &  A,
const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, typename CrsMatrixType::device_type > &  lclRowInds 
)

For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on device, using the default execution space instance.

Template Parameters
CrsMatrixTypeSpecialization of Tpetra::CrsMatrix.
Parameters
A[in/out] Sparse matrix to modify.
lclRowInds[in] Local indices of the rows of A to modify.

Definition at line 222 of file Tpetra_applyDirichletBoundaryCondition.hpp.

template<class CrsMatrixType >
void Tpetra::applyDirichletBoundaryConditionToLocalMatrixRows ( CrsMatrixType &  A,
const Kokkos::View< typename CrsMatrixType::local_ordinal_type *, Kokkos::HostSpace > &  lclRowInds 
)

For all k in [0, lclRowInds.extent(0)), set local row lclRowInds[k] of A to have 1 on the diagonal and 0 elsewhere. Run on host, using the default host execution space instance.

Template Parameters
CrsMatrixTypeSpecialization of Tpetra::CrsMatrix.
Parameters
A[in/out] Sparse matrix to modify.
lclRowInds[in] Local indices of the rows of A to modify.

Definition at line 234 of file Tpetra_applyDirichletBoundaryCondition.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::blockCrsMatrixWriter ( BlockCrsMatrix< Scalar, LO, GO, Node > const &  A,
std::string const &  fileName 
)

Helper function to write a BlockCrsMatrix. Calls the 3-argument version.

Definition at line 59 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::blockCrsMatrixWriter ( BlockCrsMatrix< Scalar, LO, GO, Node > const &  A,
std::string const &  fileName,
Teuchos::ParameterList const &  params 
)

Helper function to write a BlockCrsMatrix. Calls the 3-argument version.

Definition at line 67 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::blockCrsMatrixWriter ( BlockCrsMatrix< Scalar, LO, GO, Node > const &  A,
std::ostream &  os 
)

Helper function to write a BlockCrsMatrix. Calls the 3-argument version.

Definition at line 74 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::blockCrsMatrixWriter ( BlockCrsMatrix< Scalar, LO, GO, Node > const &  A,
std::ostream &  os,
Teuchos::ParameterList const &  params 
)

Helper function to write a BlockCrsMatrix.

Writes the block matrix to the specified ostream in point form. The following parameter list options are available:

  • "always use parallel algorithm" : on one process, this forces the use of the parallel strip-mining algorithm (default=false)
  • "print MatrixMarket header" : if false, don't print the MatrixMarket header (default=true)
  • "precision" : precision to be used in printing matrix entries (default=C++ default)
  • "zero-based indexing" : if true, print the matrix with 0-based indexing (default=false)

Definition at line 80 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::writeMatrixStrip ( BlockCrsMatrix< Scalar, LO, GO, Node > const &  A,
std::ostream &  os,
Teuchos::ParameterList const &  params 
)

Helper function called by blockCrsMatrixWriter.

This function should not be called directly.

Definition at line 198 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > Tpetra::convertToBlockCrsMatrix ( const Tpetra::CrsMatrix< Scalar, LO, GO, Node > &  pointMatrix,
const LO &  blockSize 
)

Non-member constructor that creates a BlockCrsMatrix from an existing point CrsMatrix.

This function accepts an already constructed point version of the block matrix. Assumptions:

  • All point entries in a logical block must be stored in the CrsMatrix, even if the values are zero.
  • Point rows corresponding to a particular mesh node must be stored consecutively.

Definition at line 294 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class LO , class GO , class Node >
Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > Tpetra::createMeshMap ( LO const &  blockSize,
const Tpetra::Map< LO, GO, Node > &  pointMap 
)

Helper function to generate a mesh map from a point map. Important! It's assumed that point GIDs associated with a single mesh GID appear consecutively in pointMap.

Definition at line 418 of file Tpetra_BlockCrsMatrix_Helpers_def.hpp.

template<class ViewType , class CoefficientType , class LayoutType = typename ViewType::array_layout, class IndexType = int, const int rank = ViewType::rank>
KOKKOS_INLINE_FUNCTION void Tpetra::SCAL ( const CoefficientType &  alpha,
const ViewType &  x 
)

x := alpha*x, where x is either rank 1 (a vector) or rank 2 (a matrix).

Definition at line 664 of file Tpetra_BlockView.hpp.

template<class ViewType , class InputType , class LayoutType = typename ViewType::array_layout, class IndexType = int, const int rank = ViewType::rank>
KOKKOS_INLINE_FUNCTION void Tpetra::FILL ( const ViewType &  x,
const InputType &  val 
)

Set every entry of x to val.

Definition at line 676 of file Tpetra_BlockView.hpp.

template<class CoefficientType , class ViewType1 , class ViewType2 , class LayoutType1 = typename ViewType1::array_layout, class LayoutType2 = typename ViewType2::array_layout, class IndexType = int, const int rank = ViewType1::rank>
KOKKOS_INLINE_FUNCTION void Tpetra::AXPY ( const CoefficientType &  alpha,
const ViewType1 &  x,
const ViewType2 &  y 
)

y := y + alpha * x (dense vector or matrix update)

This function follows the BLAS convention that if alpha == 0, then it does nothing. (This matters only if x contains Inf or NaN values.)

Definition at line 694 of file Tpetra_BlockView.hpp.

template<class ViewType1 , class ViewType2 , class LayoutType1 = typename ViewType1::array_layout, class LayoutType2 = typename ViewType2::array_layout, class IndexType = int, const int rank = ViewType1::rank>
KOKKOS_INLINE_FUNCTION void Tpetra::COPY ( const ViewType1 &  x,
const ViewType2 &  y 
)

Deep copy x into y, where x and y are either rank 1 (vectors) or rank 2 (matrices) with the same dimension(s).

Parameters
x[in] The input vector / matrix.
y[out] The output vector / matrix.

We put the output argument last, because that's what the BLAS functions _COPY (replace _ with "S", "D", "C", or "Z") do.

Definition at line 720 of file Tpetra_BlockView.hpp.

template<class VecType1 , class BlkType , class VecType2 , class CoeffType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void Tpetra::GEMV ( const CoeffType &  alpha,
const BlkType &  A,
const VecType1 &  x,
const VecType2 &  y 
)

y := y + alpha * A * x (dense matrix-vector multiply)

Parameters
alpha[in] Coefficient by which to multiply A*x (this does NOT necessarily follow BLAS rules; the caller is responsible for checking whether alpha == 0 and implementing BLAS rules in that case).
A[in] Small dense matrix (must have rank 2)
x[in] Small dense vector input (must have rank 1 and at least as many rows as A has columns)
y[in/out] Small dense vector output (must have rank 1 and at least as many rows as A has rows)

Definition at line 746 of file Tpetra_BlockView.hpp.

template<class ViewType1 , class ViewType2 , class ViewType3 , class CoefficientType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void Tpetra::GEMM ( const char  transA[],
const char  transB[],
const CoefficientType &  alpha,
const ViewType1 &  A,
const ViewType2 &  B,
const CoefficientType &  beta,
const ViewType3 &  C 
)

Small dense matrix-matrix multiply: C := alpha*A*B + beta*C

Template Parameters
ViewType1Type of the first matrix input A.
ViewType2Type of the second matrix input B.
ViewType3Type of the third matrix input/output C.
CoefficientTypeType of the scalar coefficients alpha and beta.
IndexTypeType of the index used in for loops; defaults to int.

Definition at line 768 of file Tpetra_BlockView.hpp.

template<class LittleBlockType , class LittleVectorType >
KOKKOS_INLINE_FUNCTION void Tpetra::GETF2 ( const LittleBlockType &  A,
const LittleVectorType &  ipiv,
int &  info 
)

Computes A = P*L*U.

Definition at line 907 of file Tpetra_BlockView.hpp.

template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void Tpetra::GETRS ( const char  mode[],
const LittleBlockType &  A,
const LittleIntVectorType &  ipiv,
const LittleScalarVectorType &  B,
int &  info 
)

Solve the linear system(s) AX=B, using the result of GETRF or GETF2.

Warning
We have not implemented transpose yet, or multiple right-hand sides.

Definition at line 1141 of file Tpetra_BlockView.hpp.

template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void Tpetra::GETRI ( const LittleBlockType &  A,
const LittleIntVectorType &  ipiv,
const LittleScalarVectorType &  work,
int &  info 
)

Compute inverse of A, using result of GETRF or GETF2.

Template Parameters
LittleBlockTypeType of dense matrix A
LittleBlockTypeType of 1-D pivot array ipiv
LittleScalarVectorTypeType of 1-D work array work
Parameters
A[in/out] On input: output matrix resulting from running GETRF or GETF2 on a square matrix A. On output: inverse of the original matrix A.
ipiv[in] Pivot array from the LU factorization.
work[out] Temporary workspace; must be at least as long as the number of rows in A.
info[out] On output, 0 if the routine was successful, else nonzero.

Definition at line 1170 of file Tpetra_BlockView.hpp.

void Tpetra::setCombineModeParameter ( Teuchos::ParameterList &  plist,
const std::string &  paramName 
)

Set CombineMode parameter in a Teuchos::ParameterList.

If you are constructing a Teuchos::ParameterList with a CombineMode parameter, set the parameter by using this function. This will use a special feature of Teuchos – custom parameter list validation – so that users can specify CombineMode values by string, rather than enum value. The strings are the same as the enum names: "ADD", "INSERT", "REPLACE", "ABSMAX", and "ZERO". They are not case sensitive.

Using this function to set a CombineMode parameter will ensure that the XML serialization of the resulting Teuchos::ParameterList will refer to the CombineMode enum values using human-readable string names, rather than raw integers.

Parameters
plist[out] Teuchos::ParameterList to which you want to add the Tpetra::CombineMode parameter.
paramName[in] String name to use for the parameter. For example, you might wish to call the parameter "Combine Mode", "Tpetra::CombineMode", or "combine mode". The parameter's name is case sensitive, even though the string values are not.

Definition at line 48 of file Tpetra_CombineMode.cpp.

std::string Tpetra::combineModeToString ( const CombineMode  combineMode)

Human-readable string representation of the given CombineMode.

Definition at line 89 of file Tpetra_CombineMode.cpp.

template<class SC , class LO , class GO , class NT >
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > Tpetra::computeRowOneNorms ( const Tpetra::RowMatrix< SC, LO, GO, NT > &  A)

Compute global row one-norms ("row sums") of the input sparse matrix A, in a way suitable for one-sided (left only) equilibration.

Note
This function is collective over A's communicator, and may need to communicate, depending on A's Maps.
Parameters
A[in] The input sparse matrix A.
Returns
Input to leftAndOrRightScaleCrsMatrix (which see). The result is only safe to use for left scaling, not for right scaling.

Definition at line 965 of file Tpetra_computeRowAndColumnOneNorms_def.hpp.

template<class SC , class LO , class GO , class NT >
Details::EquilibrationInfo< typename Kokkos::ArithTraits< SC >::val_type, typename NT::device_type > Tpetra::computeRowAndColumnOneNorms ( const Tpetra::RowMatrix< SC, LO, GO, NT > &  A,
const bool  assumeSymmetric 
)

Compute global row and column one-norms ("row sums" and "column sums") of the input sparse matrix A, in a way suitable for two-sided (left and right) equilibration.

Note
For AztecOO users: If you set assumeSymmetric=true, this function should behave like setting the AZ_scaling option to AZ_sym_row_sum.
This function is collective over A's communicator, and may need to communicate, depending on A's Maps.

For the nonsymmetric case, this function works like a sparse version of LAPACK's DGEEQU routine, except that it uses one norms (sums of absolute values) instead of infinity norms (maximum absolute value). The resulting row and column scaling is NOT symmetric. For the symmetric case, this function computes the row norms and uses those for the column norms. The resulting scaling is symmetric IF you take square roots.

Parameters
A[in] The input sparse matrix A.
assumeSymmetric[in] Whether to assume that the matrix A is (globally) symmetric. If so, don't compute row-scaled column norms separately from row norms.
Returns
Input to leftAndOrRightScaleCrsMatrix (which see).

Definition at line 979 of file Tpetra_computeRowAndColumnOneNorms_def.hpp.

bool Tpetra::isInitialized ( )

Whether Tpetra is in an initialized state.

Initialize Tpetra by calling one of the versions of initialize(). After initialize() returns, Tpetra is initialized. Once finalize() returns, Tpetra is no longer initialized.

Definition at line 195 of file Tpetra_Core.cpp.

Teuchos::RCP< const Teuchos::Comm< int > > Tpetra::getDefaultComm ( )

Get Tpetra's default communicator.

Precondition
One of the Tpetra::initialize() functions has been called.
Returns
If one of the versions of initialize() was called that takes a default communicator, this function returns that communicator. Otherwise, this function returns MPI_COMM_WORLD (wrapped in a Teuchos wrapper) if Trilinos was built with MPI enabled, or a Teuchos::SerialComm instance otherwise.

Definition at line 199 of file Tpetra_Core.cpp.

void Tpetra::initialize ( int *  argc,
char ***  argv 
)

Initialize Tpetra.

This initializes the following if they have not already been initialized:

  • MPI (if Trilinos was built with MPI enabled)
  • Kokkos

If Trilinos was built with MPI enabled, this function sets the default communicator to MPI_COMM_WORLD (wrapped in a Teuchos wrapper). Otherwise, it sets the default communicator to a Teuchos::SerialComm instance.

Parameters
argc[in/out] Same as first argument of MPI_Init()
argv[in/out] Same as second argument of MPI_Init()

The argc and argv arguments are both passed by pointer, in order to match MPI_Init's interface. MPI_Init() reserves the right to modify command-line arguments, e.g., by reading and removing those that pertain to MPI. Thus, in main(), one would write

Tpetra::initialize (&argc, &argc);

Definition at line 228 of file Tpetra_Core.cpp.

void Tpetra::initialize ( int *  argc,
char ***  argv,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Initialize Tpetra.

Warning
It is NOT legal to create a Teuchos::MpiComm until after MPI has been initialized. Thus, you may not call this function with a Teuchos::MpiComm instance unless MPI has been initialized.

This initializes the following if they have not already been initialized:

  • MPI (if Trilinos was built with MPI enabled)
  • Kokkos
Parameters
argc[in/out] Same as first argument of MPI_Init()
argv[in/out] Same as second argument of MPI_Init()
comm[in] Tpetra's default communicator, wrapped in a Teuchos wrapper. This may be either a Teuchos::MpiComm or a Teuchos::SerialComm instance.

The argc and argv arguments are both passed by pointer, in order to match MPI_Init's interface. MPI_Init() reserves the right to modify command-line arguments, e.g., by reading and removing those that pertain to MPI. Thus, in main(), one would write

Teuchos::RCP<const Teuchos::Comm<int> > comm = ...; // whatever you want
Tpetra::initialize (&argc, &argc, comm);

Definition at line 289 of file Tpetra_Core.cpp.

void Tpetra::finalize ( )

Finalize Tpetra.

If Tpetra::initialize initialized Kokkos, finalize Kokkos. If Tpetra::initialize initialized MPI, finalize MPI. Don't call this unless you first call Tpetra::initialize.

If you (the user) initialized Kokkos resp. MPI before Tpetra::initialize was called, then this function does NOT finalize Kokkos resp. MPI. In that case, you (the user) are responsible for finalizing Kokkos resp. MPI.

Definition at line 306 of file Tpetra_Core.cpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > createCrsGraph ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &  map,
size_t  maxNumEntriesPerRow = 0,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

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

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

Definition at line 2344 of file Tpetra_CrsGraph_decl.hpp.

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

Nonmember CrsGraph constructor that fuses Import and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

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

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

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

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

Definition at line 2409 of file Tpetra_CrsGraph_decl.hpp.

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

Nonmember CrsGraph constructor that fuses Import and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

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

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

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

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

Definition at line 2478 of file Tpetra_CrsGraph_decl.hpp.

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

Nonmember CrsGraph constructor that fuses Export and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

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

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

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

Definition at line 2533 of file Tpetra_CrsGraph_decl.hpp.

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

Nonmember CrsGraph constructor that fuses Export and fillComplete().

Template Parameters
CrsGraphTypeA specialization of CrsGraph.

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

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

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

Definition at line 2585 of file Tpetra_CrsGraph_decl.hpp.

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > Tpetra::importAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  importer,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Nonmember CrsMatrix constructor that fuses Import and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

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

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.

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

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

Definition at line 4489 of file Tpetra_CrsMatrix_decl.hpp.

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > Tpetra::importAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  rowImporter,
const Import< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  domainImporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)

Nonmember CrsMatrix constructor that fuses Import and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

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

Fusing redistribution and fillComplete() exposes potential optimizations. For example, it may make constructing the column Map faster, and it may avoid intermediate unoptimized storage in the destination CrsMatrix. These optimizations may improve performance for specialized kernels like sparse matrix-matrix multiply, as well as for redistributing data after doing load balancing.

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

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

Definition at line 4508 of file Tpetra_CrsMatrix_decl.hpp.

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > Tpetra::exportAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  exporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap = Teuchos::null,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap = Teuchos::null,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Nonmember CrsMatrix constructor that fuses Export and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

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

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

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

Definition at line 4530 of file Tpetra_CrsMatrix_decl.hpp.

template<class CrsMatrixType >
Teuchos::RCP< CrsMatrixType > Tpetra::exportAndFillCompleteCrsMatrix ( const Teuchos::RCP< const CrsMatrixType > &  sourceMatrix,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  rowExporter,
const Export< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > &  domainExporter,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  domainMap,
const Teuchos::RCP< const Map< typename CrsMatrixType::local_ordinal_type, typename CrsMatrixType::global_ordinal_type, typename CrsMatrixType::node_type > > &  rangeMap,
const Teuchos::RCP< Teuchos::ParameterList > &  params 
)

Nonmember CrsMatrix constructor that fuses Export and fillComplete().

Template Parameters
CrsMatrixTypeA specialization of CrsMatrix.

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

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

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

Definition at line 4549 of file Tpetra_CrsMatrix_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createCrsMatrix ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const size_t  maxNumEntriesPerRow = 0,
const Teuchos::RCP< Teuchos::ParameterList > &  params = Teuchos::null 
)

Create an empty CrsMatrix given a row map and a single integer upper bound on the number of stored entries per row.

Definition at line 4475 of file Tpetra_CrsMatrix_decl.hpp.

template<class OpScalar , class MatScalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< CrsMatrixMultiplyOp< OpScalar, MatScalar, LocalOrdinal, GlobalOrdinal, Node > > createCrsMatrixMultiplyOp ( const Teuchos::RCP< const CrsMatrix< MatScalar, LocalOrdinal, GlobalOrdinal, Node > > &  A)

Non-member function to create a CrsMatrixMultiplyOp.

The function has the same template parameters of CrsMatrixMultiplyOp.

Parameters
A[in] The CrsMatrix instance to wrap in an CrsMatrixMultiplyOp.
Returns
The CrsMatrixMultiplyOp wrapper for the given CrsMatrix.

Definition at line 1111 of file Tpetra_CrsMatrixMultiplyOp.hpp.

template<class DistObjectType >
void Tpetra::removeEmptyProcessesInPlace ( Teuchos::RCP< DistObjectType > &  input,
const Teuchos::RCP< const Map< typename DistObjectType::local_ordinal_type, typename DistObjectType::global_ordinal_type, typename DistObjectType::node_type > > &  newMap 
)

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

Template Parameters
DistObjectTypeA specialization of DistObject.
Warning
This method is ONLY for use by experts. The fact that the documentation of this method starts with a "Vocabulary" section should give you proper respect for the complicated semantics of this method in a parallel MPI run.
We make NO promises of backwards compatibility. This method may change or disappear at any time.

Vocabulary:

  • The Map returned by input->getMap() on input to this method is the "original Map."
  • The communicator returned by input->getComm() on input to this method is the "original communicator."
  • All processes in the original communicator which contain zero elements in the original Map are "excluded processes."
  • All other processes in the original communicator are "included processes."

Preconditions:

  • The nonnull object input is distributed over the original Map.
  • The input Map newMap must be the same as the result of calling removeEmptyProcesses() on the original Map.
  • On excluded processes, newMap must be Teuchos::null. (This is what getMap()->removeEmptyProcesses() returns anyway on excluded processes.)

This method has collective semantics over the original communicator. 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. The input DistObject input will be nonnull on return. On excluded processes, free any data in input that need freeing, do any other work that needs to be done to restore correct semantics, and set input to null before returning.

The two-argument version of this function is useful if you have already precomputed the new Map that excludes processes with zero elements. For example, you might want to apply this Map to several different MultiVector instances. The one-argument version of this function is useful if you want the DistObject to compute the new Map itself, because you only plan to use it for that one DistObject instance.

Here is a sample use case. Suppose that input is some subclass of DistObject, like MultiVector, CrsGraph, or CrsMatrix. Suppose also that map_type is the corresponding specialization of Map.

RCP<const map_type> origRowMap = input->getMap ();
RCP<const map_type> newRowMap = origRowMap->removeEmptyProcesses ();
removeEmptyProcessesInPlace (input, newRowMap);
// Either (both the new Map and input are null), or
// (both the new Map and input are not null).
assert ((newRowMap.is_null () && input.is_null ()) ||
(! newRowMap.is_null () && ! input.is_null ()));

Warning
On excluded processes, calling this function invalidates any other references to the input DistObject input. Calling any methods (other than the destructor) on the input on excluded processes has undefined behavior in that case, and may result in deadlock.
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.
To implementers of DistObject subclasses: On exit, the only method of this object which is safe to call on excluded processes is the destructor, or this method with the original Map. This implies that subclasses' destructors must not contain communication operations.

template<class DistObjectType >
void Tpetra::removeEmptyProcessesInPlace ( Teuchos::RCP< DistObjectType > &  input)

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

Template Parameters
DistObjectTypeA specialization of DistObject.
Warning
This method is ONLY for use by experts.
We make NO promises of backwards compatibility. This method may change or disappear at any time.

This method behaves just like the two-argument version of removeEmptyProcessesInPlace(), except that it first calls removeEmptyProcesses() on the input DistObject's Map to compute the new Map.

The two-argument version of this function is useful if you have already precomputed the new Map that excludes processes with zero elements. For example, you might want to apply this Map to several different MultiVector instances. The one-argument version of this function is useful if you want the DistObject to compute the new Map itself, because you only plan to use it for that one DistObject instance.

Here is a sample use case. Suppose that input is some subclass of DistObject, like MultiVector, CrsGraph, or CrsMatrix. Suppose also that map_type is the corresponding specialization of Map.

RCP<const map_type> newRowMap;
if (! input.is_null ()) {
newRowMap = input->getMap ();
}
// Either (both the new Map and input are null), or
// (both the new Map and input are not null).
assert ((newRowMap.is_null () && input.is_null ()) ||
(! newRowMap.is_null () && ! input.is_null ()));

Definition at line 1307 of file Tpetra_DistObject_def.hpp.

Teuchos::Array< std::string > Tpetra::distributorSendTypes ( )

Valid values for Distributor's "Send type" parameter.

This is mainly useful as an implementation detail of Distributor. You may use it if you would like a programmatic way to get all possible values of the "Send type" parameter of Distributor.

Definition at line 94 of file Tpetra_Distributor.cpp.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > createExport ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  src,
const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  tgt 
)

Nonmember "constructor" for Export objects.

Create an Export object from the given source and target Maps.

Precondition
src != null
tgt != null
Returns
The Export object. If src == tgt, returns null.

(Debug mode: throws std::runtime_error if one of src or tgt is null.)

Definition at line 290 of file Tpetra_Export_decl.hpp.

template<class ExecutionSpace , class GlobalDataStructure , class UserFunctionType >
void Tpetra::for_each ( const char  kernelLabel[],
ExecutionSpace  execSpace,
GlobalDataStructure &  X,
UserFunctionType  f 
)

Apply a function entrywise to each local entry of a Tpetra global data structure, analogously to std::for_each.

Generically, the function f takes the current entry by nonconst reference, in read-write fashion. Overloads for specific Tpetra data structures may permit functions that take other arguments.

Parameters
kernelLabel[in] Kernel label for Kokkos Profiling.
execSpace[in] Kokkos execution space on which to run.
X[in/out] Tpetra global data structure to modify.
f[in] Function to apply to each entry of X.

There is no overload without a kernel label. You WILL label your kernels. If you bother to use Kokkos, that means you care about performance, so we need to be able to measure it.

To implement this for new GlobalDataStructure types, specialize Details::ForEach (see below).

Definition at line 125 of file Tpetra_for_each.hpp.

template<class GlobalDataStructure , class UserFunctionType >
void Tpetra::for_each ( const char  kernelLabel[],
GlobalDataStructure &  X,
UserFunctionType  f 
)

Overload of for_each (see above) that runs on X's default Kokkos execution space.

Parameters
kernelLabel[in] Kernel label for Kokkos Profiling.
X[in/out] Tpetra global data structure to modify.
f[in] Function to apply entrywise to X (could have different signatures; see above).

There is no overload without a kernel label. You WILL label your kernels. If you bother to use Kokkos, that means you care about performance, so we need to be able to measure it.

To implement this for new GlobalDataStructure types, specialize Details::ForEach (see below).

Definition at line 141 of file Tpetra_for_each.hpp.

template<class SC , class LO , class GO , class NT >
std::shared_ptr< ::Tpetra::Details::CommRequest> Tpetra::idot ( typename::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *  resultRaw,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  X,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  Y 
)

Nonblocking dot product, with either Tpetra::MultiVector or Tpetra::Vector inputs, and raw pointer or raw array output.

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)).
Parameters
resultRaw[out] Output; raw pointer or raw array to the return value(s). It is only valid to read this after calling wait() on the return value. It must be legal to write to the first std::max(X.getNumVectors(), Y.getNumVectors()) entries of this array. This memory must be accessible from the host CPU.
X[in] First input Tpetra::MultiVector or Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as Y. If this is a Tpetra::MultiVector, then this must either have the same number of columns as Y, or have one column.
Y[in] Second input Tpetra::MultiVector or Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as X. If this is a Tpetra::MultiVector, then this must either have the same number of columns as X, or have one column.
Returns
Pointer to an object representing the nonblocking collective (communication operation). Call wait() on this object to complete the collective. After calling wait(), you may read the result.
Template Parameters
SCSame as the first template parameter of Tpetra::Vector.
LOSame as the second template parameter of Tpetra::Vector.
GOSame as the third template parameter of Tpetra::Vector.
NTSame as the fourth template parameter of Tpetra::Vector.

In this version of the function, the dot product result goes into an array (just a raw pointer). The dot_type typedef is the type of a dot product result, for a Tpetra::Vector whose entries have type SC (the "Scalar" type). For most SC types, dot_type == SC. However, once you start dipping into more interesting Scalar types, such as those found in the Sacado or Stokhos packages, dot_type may be a different type. Most users should not have to worry about this, but Tpetra developers may need to worry about this.

Definition at line 373 of file Tpetra_idot.hpp.

template<class SC , class LO , class GO , class NT >
std::shared_ptr< ::Tpetra::Details::CommRequest> Tpetra::idot ( const Kokkos::View< typename::Tpetra::Vector< SC, LO, GO, NT >::dot_type, typename::Tpetra::Vector< SC, LO, GO, NT >::device_type > &  result,
const ::Tpetra::Vector< SC, LO, GO, NT > &  X,
const ::Tpetra::Vector< SC, LO, GO, NT > &  Y 
)

Nonblocking dot product, with Tpetra::Vector inputs, and rank-0 (single value) Kokkos::View output.

This function computes result() = dot(X,Y).

Parameters
result[out] Output; rank-0 Kokkos::View of the return value. It is only valid to read this after calling wait() on the return value.
X[in] First input Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as Y.
Y[in] Second input Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as X.
Returns
Pointer to an object representing the nonblocking collective (communication operation). Call wait() on this object to complete the collective. After calling wait(), you may read the result.
Template Parameters
SCSame as the first template parameter of Tpetra::Vector.
LOSame as the second template parameter of Tpetra::Vector.
GOSame as the third template parameter of Tpetra::Vector.
NTSame as the fourth template parameter of Tpetra::Vector.

In this version of the function, the dot product result goes into a rank-0 ("zero-dimensional") Kokkos::View. Rank-0 Views just view a single value. We prefer that you use the versions of idot() that take a Kokkos::View as their output argument, because this ensures that the output will still exist (not be deallocated or fall out of scope). The versions of idot() that take a raw pointer cannot promise that the memory will continue to exist until the dot product is done.

The dot_type typedef is the type of a dot product result, for a Tpetra::Vector whose entries have type SC (the "Scalar" type). For most SC types, dot_type == SC. However, once you start dipping into more interesting Scalar types, such as those found in the Sacado or Stokhos packages, dot_type may be a different type. Most users should not have to worry about this, but Tpetra developers may need to worry about this.

Definition at line 480 of file Tpetra_idot.hpp.

template<class SC , class LO , class GO , class NT >
std::shared_ptr< ::Tpetra::Details::CommRequest> Tpetra::idot ( const Kokkos::View< typename::Tpetra::MultiVector< SC, LO, GO, NT >::dot_type *, typename::Tpetra::MultiVector< SC, LO, GO, NT >::device_type > &  result,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  X,
const ::Tpetra::MultiVector< SC, LO, GO, NT > &  Y 
)

Nonblocking dot product, with Tpetra::MultiVector inputs, and rank-1 (one-dimensional array) Kokkos::View output.

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, result[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, result[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, result[j] = dot(X(:,j), Y(:,0)).
Parameters
result[out] Output; rank-1 Kokkos::View. It is only valid to read the entries of this after calling wait() on the return value. This must have length std::max(X.getNumVectors(), Y.getNumVectors()).
X[in] First input Tpetra::MultiVector or Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as Y. If this is a Tpetra::MultiVector, then this must either have the same number of columns as Y, or have one column.
Y[in] Second input Tpetra::MultiVector or Tpetra::Vector. This must have same number of rows (globally, and on each (MPI) process) as X. If this is a Tpetra::MultiVector, then this must either have the same number of columns as X, or have one column.
Returns
Pointer to an object representing the nonblocking collective (communication operation). Call wait() on this object to complete the collective. After calling wait(), you may read the results.
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.

Compute the dot product of each column of X, with each corresponding column of Y. If X has a single column, then compute the dot product of X with each column of Y in turn. If Y has a single column, then compute the dot product of each column of X in turn with Y.

In this version of the function, the dot product results go into a rank-1 (one-dimensional) Kokkos::View. We prefer that you use the versions of idot() that take a Kokkos::View as their output argument, because this ensures that the output will still exist (not be deallocated or fall out of scope). The versions of idot() that take a raw pointer cannot promise that the memory will continue to exist until the dot product is done.

The dot_type typedef is the type of a dot product result, for a Tpetra::Vector whose entries have type SC (the "Scalar" type). For most SC types, dot_type == SC. However, once you start dipping into more interesting Scalar types, such as those found in the Sacado or Stokhos packages, dot_type may be a different type. Most users should not have to worry about this, but Tpetra developers may need to worry about this.

Definition at line 578 of file Tpetra_idot.hpp.

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

Nonmember constructor for Import.

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

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

Definition at line 537 of file Tpetra_Import_decl.hpp.

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

Nonmember constructor for Import that takes a ParameterList.

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

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

Definition at line 564 of file Tpetra_Import_decl.hpp.

template<class SC , class LO , class GO , class NT >
void Tpetra::leftAndOrRightScaleCrsMatrix ( Tpetra::CrsMatrix< SC, LO, GO, NT > &  A,
const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &  rowScalingFactors,
const Kokkos::View< const typename Kokkos::ArithTraits< SC >::mag_type *, typename NT::device_type > &  colScalingFactors,
const bool  leftScale,
const bool  rightScale,
const bool  assumeSymmetric,
const EScaling  scaling 
)

Left-scale and/or right-scale (in that order) the entries of the input Tpetra::CrsMatrix A.

Parameters
A[in/out] The sparse matrix A to scale. It must have a valid KokkosSparse::CrsMatrix. This is true if fillComplete has been called on it at least once, or if the matrix was created with a local sparse matrix.
rowScalingFactors[in] Use Details::EquilibrationInfo::rowNorms.
colScalingFactors[in] If assumeSymmetric, use Details::EquilibrationInfo::colNorms, else use Details::EquilibrationInfo::rowScaledColNorms.
leftScale[in] Whether to left-scale A. Left scaling happens first.
rightScale[in] Whether to right-scale A. Right scaling happens last.
Whetherto assume symmetric scaling, that is, whether to take square roots of scaling factors before scaling.
scaling[in] If SCALING_DIVIDE, "scale" means "divide by"; if SCALING_MULTIPLY, it means "multiply by."

Definition at line 63 of file Tpetra_leftAndOrRightScaleCrsMatrix_def.hpp.

template<class SC , class LO , class GO , class NT >
void Tpetra::leftAndOrRightScaleCrsMatrix ( Tpetra::CrsMatrix< SC, LO, GO, NT > &  A,
const Tpetra::Vector< typename Kokkos::ArithTraits< SC >::mag_type, LO, GO, NT > &  rowScalingFactors,
const Tpetra::Vector< typename Kokkos::ArithTraits< SC >::mag_type, LO, GO, NT > &  colScalingFactors,
const bool  leftScale,
const bool  rightScale,
const bool  assumeSymmetric,
const EScaling  scaling 
)

Left-scale and/or right-scale (in that order) the entries of the input Tpetra::CrsMatrix A.

Parameters
A[in/out] The sparse matrix A to scale. It must have a valid KokkosSparse::CrsMatrix. This is true if fillComplete has been called on it at least once, or if the matrix was created with a local sparse matrix.
rowScalingFactors[in] The row scaling factors; must be distributed according to the matrix's row Map. This function does NOT redistribute the Vector for you.
colScalingFactors[in] The column scaling factors; must be distributed according to the matrix's column Map. This function does NOT redistribute the Vector for you.
leftScale[in] Whether to left-scale A. Left scaling happens first.
rightScale[in] Whether to right-scale A. Right scaling happens last.
Whetherto assume symmetric scaling, that is, whether to take square roots of scaling factors before scaling.
scaling[in] If SCALING_DIVIDE, "scale" means "divide by"; if SCALING_MULTIPLY, it means "multiply by."

Definition at line 124 of file Tpetra_leftAndOrRightScaleCrsMatrix_def.hpp.

template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createLocalMap ( const size_t  numElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Nonmember constructor for a locally replicated Map with the default Kokkos Node.

This method returns a Map instantiated on the default Kokkos Node type. The Map is configured to use zero-based indexing.

Parameters
numElements[in] Number of elements on each process. Each process gets the same set of elements, namely 0, 1, ..., numElements - 1.
comm[in] The Map's communicator.
template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createLocalMapWithNode ( const size_t  numElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Nonmember constructor for a locally replicated Map with a specified Kokkos Node.

This method returns a Map instantiated on the given Kokkos Node instance. The Map is configured to use zero-based indexing.

Parameters
numElements[in] Number of elements on each process. Each process gets the same set of elements, namely 0, 1, ..., numElements - 1.
comm[in] The Map's communicator.
template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createUniformContigMap ( const global_size_t  numElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Non-member constructor for a uniformly distributed, contiguous Map with the default Kokkos Node.

This method returns a Map instantiated on the Kokkos default Node type. The resulting Map uses zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createUniformContigMapWithNode ( const global_size_t  numElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Non-member constructor for a uniformly distributed, contiguous Map with a user-specified Kokkos Node.

The resulting Map uses zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createContigMap ( const global_size_t  numElements,
const size_t  localNumElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Non-member constructor for a (potentially) non-uniformly distributed, contiguous Map using the default Kokkos::Device.

The Map is configured to use zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createContigMapWithNode ( const global_size_t  numElements,
const size_t  localNumElements,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Nonmember constructor for a (potentially) nonuniformly distributed, contiguous Map for a user-specified, possibly nondefault Kokkos Node type.

If Node is the default, use createContigMap instead. The Map is configured to use zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createNonContigMap ( const Teuchos::ArrayView< const GlobalOrdinal > &  elementList,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Nonmember constructor for a non-contiguous Map using the default Kokkos::Device type.

The Map is configured to use zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createNonContigMapWithNode ( const Teuchos::ArrayView< const GlobalOrdinal > &  elementList,
const Teuchos::RCP< const Teuchos::Comm< int > > &  comm 
)

Nonmember constructor for a noncontiguous Map with a user-specified, possibly nondefault Kokkos Node type.

If Node is the default, use createNonContigMap instead. The Map is configured to use zero-based indexing.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createOneToOne ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  M)

Nonmember constructor for a contiguous Map with user-defined weights and a user-specified, possibly nondefault Kokkos Node type.

The Map is configured to use zero-based indexing.Creates a one-to-one version of the given Map where each GID lives on only one process.

template<class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createOneToOne ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  M,
const ::Tpetra::Details::TieBreak< LocalOrdinal, GlobalOrdinal > &  tie_break 
)

Creates a one-to-one version of the given Map where each GID lives on only one process. The given TieBreak object specifies the rule to break ties.

template<class DS , class DL , class DG , class DN , class SS , class SL , class SG , class SN >
void Tpetra::deep_copy ( MultiVector< DS, DL, DG, DN > &  dst,
const MultiVector< SS, SL, SG, SN > &  src 
)

Copy the contents of the MultiVector src into dst.

Precondition
The two inputs must have the same communicator.
The Map of src must be compatible with the Map of dst.
The two inputs must have the same number of columns.

Copy the contents of the MultiVector src into the MultiVector dst. ("Copy the contents" means the same thing as "deep copy.") The two MultiVectors need not necessarily have the same template parameters, but the assignment of their entries must make sense. Furthermore, their Maps must be compatible, that is, the MultiVectors' local dimensions must be the same on all processes.

This method must always be called as a collective operation on all processes over which the multivector is distributed. This is because the method reserves the right to check for compatibility of the two Maps, at least in debug mode, and throw if they are not compatible.

Definition at line 2393 of file Tpetra_MultiVector_decl.hpp.

template<class ST , class LO , class GO , class NT >
MultiVector< ST, LO, GO, NT > Tpetra::createCopy ( const MultiVector< ST, LO, GO, NT > &  src)

Return a deep copy of the given MultiVector.

Note
MultiVector's constructor returns a shallow copy of its input, by default. If you want a deep copy, use the two-argument copy constructor with Teuchos::Copy as the second argument, or call this function (createCopy).

Definition at line 4679 of file Tpetra_MultiVector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createMultiVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map,
const size_t  numVectors 
)

Nonmember MultiVector "constructor": Create a MultiVector from a given Map.

Parameters
map[in] Map describing the distribution of rows of the resulting MultiVector.
numVectors[in] Number of columns of the resulting MultiVector.
template<class ST , class LO , class GO , class NT >
void Tpetra::deep_copy ( MultiVector< ST, LO, GO, NT > &  dst,
const MultiVector< ST, LO, GO, NT > &  src 
)

Specialization of deep_copy for MultiVector objects with the same template parameters.

Definition at line 2380 of file Tpetra_MultiVector_decl.hpp.

template<class SC , class LO , class GO , class NT >
LO Tpetra::replaceDiagonalCrsMatrix ( ::Tpetra::CrsMatrix< SC, LO, GO, NT > &  matrix,
const ::Tpetra::Vector< SC, LO, GO, NT > &  newDiag 
)

Replace diagonal entries of an input Tpetra::CrsMatrix matrix with values given in newDiag.

Template Parameters
ScalarThe type of the numerical entries of the matrix. (You can use real-valued or complex-valued types here, unlike in Epetra, where the scalar type is always double.)
LocalOrdinalThe type of local indices. See the documentation of Map for requirements.
GlobalOrdinalThe type of global indices. See the documentation of Map for requirements.
NodeThe Kokkos Node type. See the documentation of Map for requirements.
Parameters
in/out]matrix Tpetra::CrsMatrix to be modified [in] newDiag Tpetra::Vector with new values for the diagonal
Returns
Local number of successfully replaced diagonal entries
template<class ExecutionSpace , class GlobalDataStructure , class UnaryFunctionType >
void Tpetra::transform ( const char  kernelLabel[],
ExecutionSpace  execSpace,
GlobalDataStructure &  input,
GlobalDataStructure &  output,
UnaryFunctionType  f 
)

For each local entry input_i of input, assign f(input_i) to the corresponding local entry output_i of output, analogously to std::transform.

Unary transform on the given execution space instance execSpace: output_i = f(input_i).

Parameters
kernelLabel[in] Kernel label for Kokkos Profiling.
execSpace[in] Kokkos execution space on which to run.
input[in] Tpetra global data structure to read.
output[out] Tpetra global data structure to write.
f[in] Entrywise function that takes an input object's entry by const reference, and returns a value assignable to an entry of the output object.

There is no overload without a kernel label. You WILL label your kernels. If you bother to use Kokkos, that means you care about performance, so we need to be able to measure it.

To implement this for new GlobalDataStructure types, specialize Details::Transform (see below).

Definition at line 142 of file Tpetra_transform.hpp.

template<class GlobalDataStructure , class UnaryFunctionType >
void Tpetra::transform ( const char  kernelLabel[],
GlobalDataStructure &  input,
GlobalDataStructure &  output,
UnaryFunctionType  f 
)

Overload of transform (see above) that runs on the output object's default Kokkos execution space.

Unary transform on output's default execution space: output_i = f(input_i).

Parameters
kernelLabel[in] Kernel label for Kokkos Profiling.
input[in] Tpetra global data structure to read.
output[out] Tpetra global data structure to write.
f[in] Entrywise function that takes an input object's entry by const reference, and returns a value assignable to an entry of the output object.

There is no overload without a kernel label. You WILL label your kernels. If you bother to use Kokkos, that means you care about performance, so we need to be able to measure it.

To implement this for new GlobalDataStructure types, specialize Details::ForEach (see below).

Definition at line 161 of file Tpetra_transform.hpp.

template<class ExecutionSpace , class GlobalDataStructure , class BinaryFunctionType >
void Tpetra::transform ( const char  kernelLabel[],
ExecutionSpace  execSpace,
GlobalDataStructure &  input1,
GlobalDataStructure &  input2,
GlobalDataStructure &  output,
BinaryFunctionType  f 
)

Binary transform on the given execution space execSpace: output_i = f(input1_i, input2_i).

Definition at line 183 of file Tpetra_transform.hpp.

template<typename MapType , typename KeyArgType , typename ValueArgType >
MapType::iterator Tpetra::efficientAddOrUpdate ( MapType &  m,
const KeyArgType &  k,
const ValueArgType &  v 
)

Efficiently insert or replace an entry in an std::map.

Template Parameters
MapTypeSpecialization of std::map
KeyArgTypeType of keys of the std::map
ValueArgTypeType of values of the std::map

This function is taken from Scott Meyers' "Effective STL", Item

  1. If the given std::map m already contains an entry with key k, replace its value with the given value v. Otherwise, insert (k,v) into the std::map. In both cases, return an iterator that points to the inserted or updated entry.

Definition at line 236 of file Tpetra_Util.hpp.

template<class IT1 , class IT2 >
void Tpetra::sort2 ( const IT1 &  first1,
const IT1 &  last1,
const IT2 &  first2 
)

Sort the first array, and apply the resulting permutation to the second array.

Sort the values in the first array (represented by the exclusive iterator range first1,last1) in ascending order. Apply the permutation resulting from the sort to the second array (represented by a starting iterator first2).

Parameters
first1A random access iterator pointing to the beginning of the first array.
last1A random access iterator pointing to the end (exclusive) of the first array.
first2A random access iterator pointing to the beginning of the second array. The second array must have no fewer elements than the first array. If the first array has N elements, then the permutation will only be applied to the first N elements of the second array.

Definition at line 532 of file Tpetra_Util.hpp.

template<class IT1 , class IT2 , class IT3 >
void Tpetra::sort3 ( const IT1 &  first1,
const IT1 &  last1,
const IT2 &  first2,
const IT3 &  first3 
)

Sort the first array, and apply the same permutation to the second and third arrays.

Sort the values in the first array (represented by the exclusive iterator range first1,last1) in ascending order. Apply the permutation resulting from the sort to the second array (represented by a starting iterator first2) and third array (represented by a starting iterator first3).

Parameters
first1A random access iterator pointing to the beginning of the first array.
last1A random access iterator pointing to the end (exclusive) of the first array.
first2A random access iterator pointing to the beginning of the second array.
first3A random access iterator pointing to the beginning of the third array.

Definition at line 566 of file Tpetra_Util.hpp.

template<class IT1 , class IT2 >
void Tpetra::merge2 ( IT1 &  indResultOut,
IT2 &  valResultOut,
IT1  indBeg,
IT1  indEnd,
IT2  valBeg,
IT2   
)

Merge values in place, additively, with the same index.

Template Parameters
IT1Iterator type for the range of indices
IT2Iterator type for the range of values

indBeg, indEnd defines a half-exclusive (does not include the end) range of indices, and valBeg, valEnd its corresponding range of values. The range of values must have the same number of entries as the range of indices. In every nondecreasing subsequence of indices, this method will merge values that have the same index, by adding the values together. When done, it assigns the new end (exclusive) of the index range to indResultOut, and the new end (exclusive) of the value range to valResultOut. (It is legal for the index range not to be sorted, but then only nondecreasing subsequences will get merged.)

For example, if the indices on input are {0, 1, 1, 3, -1, -1, -1, 0}, and their corresponding values on input are {42.0, -4.0, -3.0, 1.5, 1.0, 2.0, 3.0}, then on exit from this function, the indices are {0, 1, 3, -1, 0}, and the values are {42.0, -7.0, 1.5, 6.0, 100.0}.

On entry to the function, indResultOut may alias indEnd, and valResultOut may alias valEnd. For example, the following code is legal:

std::vector<int> ind (...);
std::vector<double> val (...);
// ... fill ind and val ...
std::vector<int>::iterator indEnd = ind.end ();
std::vector<int>::iterator valEnd = val.end ();
merge2 (indEnd, valEnd, ind.begin (), indEnd, val.begin (), valEnd);

However, the following code is not legal, because the return value of std::vector::end() cannot be modified:

std::vector<int> ind (...);
std::vector<double> val (...);
// ... fill ind and val ...
merge2 (ind.end (), val.end (), ind.begin (), ind.end (),
val.begin (), val.end ());

Definition at line 633 of file Tpetra_Util.hpp.

template<class IT1 , class IT2 , class BinaryFunction >
void Tpetra::merge2 ( IT1 &  indResultOut,
IT2 &  valResultOut,
IT1  indBeg,
IT1  indEnd,
IT2  valBeg,
IT2  valEnd,
BinaryFunction  f 
)

Merge values in place with the same index, using any associative binary function.

Template Parameters
IT1Iterator type for the range of indices
IT2Iterator type for the range of values
BinaryFunctionThe type of a function that takes two values and returns another value.

indBeg, indEnd defines a half-exclusive (does not include the end) range of indices, and valBeg, valEnd its corresponding range of values. The range of values must have the same number of entries as the range of indices. In every nondecreasing subsequence of indices, this method will merge values that have the same index, by using the given binary function. When done, it assigns the new end (exclusive) of the index range to indResultOut, and the new end (exclusive) of the value range to valResultOut. (It is legal for the index range not to be sorted, but then only nondecreasing subsequences will get merged.)

For example, if the indices on input are {0, 1, 1, 3, -1, -1, -1, 0}, their corresponding values on input are {42.0, -4.0, -3.0, 1.5, 1.0, 2.0, 3.0}, and the binary function is an instance of std::plus<double>, then on exit from this function, the indices are {0, 1, 3, -1, 0}, and the values are {42.0, -7.0, 1.5, 6.0, 100.0}.

On entry to the function, indResultOut may alias indEnd, and valResultOut may alias valEnd. For example, the following code is legal:

std::vector<int> ind (...);
std::vector<double> val (...);
// ... fill ind and val ...
std::vector<int>::iterator indEnd = ind.end ();
std::vector<int>::iterator valEnd = val.end ();
merge2 (indEnd, valEnd, ind.begin (), indEnd,
val.begin (), valEnd, std::plus<double> ());

However, the following code is not legal, because the return value of std::vector::end() cannot be modified:

std::vector<int> ind (...);
std::vector<double> val (...);
// ... fill ind and val ...
merge2 (ind.end (), val.end (), ind.begin (), ind.end (),
val.begin (), val.end (), std::plus<double> ());

Definition at line 722 of file Tpetra_Util.hpp.

template<class KeyInputIterType , class ValueInputIterType , class KeyOutputIterType , class ValueOutputIterType , class BinaryFunction >
void Tpetra::keyValueMerge ( KeyInputIterType  keyBeg1,
KeyInputIterType  keyEnd1,
ValueInputIterType  valBeg1,
ValueInputIterType  valEnd1,
KeyInputIterType  keyBeg2,
KeyInputIterType  keyEnd2,
ValueInputIterType  valBeg2,
ValueInputIterType  valEnd2,
KeyOutputIterType  keyOut,
ValueOutputIterType  valOut,
BinaryFunction  f 
)

Merge two sorted (by keys) sequences of unique (key,value) pairs by combining pairs with equal keys.

Parameters
keyBeg1[in] Start of first sequence of keys.
keyEnd1[in] End (exclusive) of first sequence of keys.
valBeg1[in] Start of first sequence of values.
valEnd1[in] End (exclusive) of first sequence of values.
keyBeg2[in] Start of second sequence of keys.
keyEnd2[in] End (exclusive) of second sequence of keys.
valBeg2[in] Start of second sequence of values.
valEnd2[in] End (exclusive) of second sequence of values.
keyOut[in/out] Output sequence of keys.
valOut[in/out] Output sequence of values.
f[in] Binary associative function to use to combine values whose keys are equal. For example, for simple replacement, use a function like std::project1st (in the SGI extensions to the STL) with both template parameters equal to the value type. For addition, use std::plus with template parameter equal to the value type.
Returns
Number of (key,value) pairs in the merged sequence.
Warning
For now, this function requires that the two input sequences be made unique by key. Later, we plan to relax that requirement.

Definition at line 790 of file Tpetra_Util.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra::createCopy ( const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &  src)

Return a deep copy of the given Vector.

Definition at line 268 of file Tpetra_Vector_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
Teuchos::RCP< Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > createVector ( const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &  map)

Nonmember Vector "constructor": Create a Vector from a given Map.

Parameters
map[in] Map describing the distribution of rows of the resulting Vector.

Definition at line 419 of file Tpetra_Vector_decl.hpp.

template<class GlobalObjectType >
Details::LocalAccess< GlobalObjectType, Details::read_only > Tpetra::readOnly ( GlobalObjectType &  G)

Declare that you want to access the given global object's local data in read-only mode.

Definition at line 685 of file Tpetra_withLocalAccess.hpp.

template<class GlobalObjectType >
Details::LocalAccess< GlobalObjectType, Details::read_only > Tpetra::readOnly ( const GlobalObjectType &  G)

Declare that you want to access the given global object's local data in read-only mode (overload for const GlobalObjectType).

Definition at line 694 of file Tpetra_withLocalAccess.hpp.

template<class GlobalObjectType >
Details::LocalAccess< GlobalObjectType, Details::write_only > Tpetra::writeOnly ( GlobalObjectType &  G)

Declare that you want to access the given global object's local data in write-only mode.

Definition at line 704 of file Tpetra_withLocalAccess.hpp.

template<class GlobalObjectType >
Details::LocalAccess< GlobalObjectType, Details::read_write > Tpetra::readWrite ( GlobalObjectType &  G)

Declare that you want to access the given global object's local data in read-and-write mode.

Definition at line 713 of file Tpetra_withLocalAccess.hpp.

template<class... LocalAccessTypes>
void Tpetra::withLocalAccess ( typename Details::ArgsToFunction< LocalAccessTypes...>::type  userFunction,
LocalAccessTypes...  localAccesses 
)

Get access to a Tpetra global object's local data.

Parameters
userFunction[in] Lambda, closure, or function that takes the local data structure by const reference.
localAccesses[in] Zero or more specializations of Details::LocalAccess. You do not create these directly; they result from applying the readOnly, writeOnly, or readWrite functions (see above) to Tpetra global objects.

A "Tpetra global object" is a Tpetra::DistObject subclass representing a data structure distributed over one or more MPI processes. Examples include Tpetra::MultiVector, Tpetra::Vector, Tpetra::CrsGraph, and Tpetra::CrsMatrix.

The type of the local data structure depends on the type of the Tpetra global object. For example, Tpetra::MultiVector's local data structure is a rank-2 unmanaged (nonowning) Kokkos::View, and Tpetra::Vector's local data structure is a rank-2 unmanaged Kokkos::View.

You may use the with_local_access_function_argument_type alias (see above) to figure out the type of the local data structure, as a function of the type of the corresponding "local access" argument. Alternately, if your function is a lambda and you have C++14, you may use const auto& ("generic lambdas").

Note
I would have preferred the LocalAccess arguments to go first. However, C++ needs the user function has to go first, else C++ can't deduce the template parameters. See the en.cppreference.com article on "parameter pack."

Definition at line 929 of file Tpetra_withLocalAccess.hpp.