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

Namespace for new Tpetra features that are not ready for public release, but are ready for evaluation by friendly expert users. More...

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

Functions

template<class Scalar , class LO , class GO , class Node >
void BlockCrsMatrix< Scalar, LO, GO, Node >copyAndPermute (const ::Tpetra::SrcDistObject &source, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs)
 
template<class Scalar , class LO , class GO , class Node >
void BlockCrsMatrix< Scalar, LO, GO, Node >packAndPrepare (const ::Tpetra::SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets, Distributor &)
 
template<class Scalar , class LO , class GO , class Node >
void BlockCrsMatrix< Scalar, LO, GO, Node >unpackAndCombine (const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t, Distributor &, const CombineMode combineMode)
 
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...
 

Detailed Description

Namespace for new Tpetra features that are not ready for public release, but are ready for evaluation by friendly expert users.

Warning
Expect header files, classes, functions, and other interfaces to change or disappear. Anything in this namespace is under active development and evaluation. Documentation may be sparse or not exist yet. Generally, unit tests will exist, but coverage may be lacking. If you understand these caveats and accept them, please feel free to take a look inside and try things out.

Function Documentation

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >copyAndPermute ( const ::Tpetra::SrcDistObject source,
const size_t  numSameIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteToLIDs,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  permuteFromLIDs 
)

Check input valid

Definition at line 2187 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >packAndPrepare ( const ::Tpetra::SrcDistObject source,
const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  exportLIDs,
Kokkos::DualView< packet_type *, buffer_device_type > &  exports,
Kokkos::DualView< size_t *, buffer_device_type >  numPacketsPerLID,
size_t &  constantNumPackets,
Distributor &   
)

Check input valid

Definition at line 2850 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::BlockCrsMatrix< Scalar, LO, GO, Node >unpackAndCombine ( const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &  importLIDs,
Kokkos::DualView< packet_type *, buffer_device_type >  imports,
Kokkos::DualView< size_t *, buffer_device_type >  numPacketsPerLID,
const size_t  ,
Distributor &  ,
const CombineMode  combineMode 
)

Check input valid

Definition at line 3103 of file Tpetra_Experimental_BlockCrsMatrix_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::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 62 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::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 70 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::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 77 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::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 83 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
void Tpetra::Experimental::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 201 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class Scalar , class LO , class GO , class Node >
Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > Tpetra::Experimental::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 297 of file Tpetra_Experimental_BlockCrsMatrix_Helpers_def.hpp.

template<class LO , class GO , class Node >
Teuchos::RCP< const Tpetra::Map< LO, GO, Node > > Tpetra::Experimental::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 421 of file Tpetra_Experimental_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::Experimental::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 672 of file Tpetra_Experimental_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::Experimental::FILL ( const ViewType &  x,
const InputType &  val 
)

Set every entry of x to val.

Definition at line 684 of file Tpetra_Experimental_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::Experimental::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 702 of file Tpetra_Experimental_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::Experimental::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 728 of file Tpetra_Experimental_BlockView.hpp.

template<class VecType1 , class BlkType , class VecType2 , class CoeffType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void Tpetra::Experimental::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 754 of file Tpetra_Experimental_BlockView.hpp.

template<class ViewType1 , class ViewType2 , class ViewType3 , class CoefficientType , class IndexType = int>
KOKKOS_INLINE_FUNCTION void Tpetra::Experimental::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 776 of file Tpetra_Experimental_BlockView.hpp.

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

Computes A = P*L*U.

Definition at line 915 of file Tpetra_Experimental_BlockView.hpp.

template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void Tpetra::Experimental::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 1149 of file Tpetra_Experimental_BlockView.hpp.

template<class LittleBlockType , class LittleIntVectorType , class LittleScalarVectorType >
KOKKOS_INLINE_FUNCTION void Tpetra::Experimental::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 1178 of file Tpetra_Experimental_BlockView.hpp.