10#ifndef TPETRA_BLOCKCRSMATRIX_DECL_HPP 
   11#define TPETRA_BLOCKCRSMATRIX_DECL_HPP 
   16#include "Tpetra_CrsGraph.hpp" 
   17#include "Tpetra_RowMatrix.hpp" 
   18#include "Tpetra_BlockMultiVector_decl.hpp" 
   21#include "KokkosSparse_BsrMatrix.hpp" 
   23#include "Tpetra_Details_MatrixApplyHelper.hpp" 
   27template <
class BlockCrsMatrixType>
 
   28Teuchos::RCP<BlockCrsMatrixType>
 
   29importAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& sourceMatrix,
 
   30                                    const Import<
typename BlockCrsMatrixType::local_ordinal_type,
 
   31                                                 typename BlockCrsMatrixType::global_ordinal_type,
 
   32                                                 typename BlockCrsMatrixType::node_type>& importer);
 
   33template <
class BlockCrsMatrixType>
 
   34Teuchos::RCP<BlockCrsMatrixType>
 
   35exportAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& sourceMatrix,
 
   36                                    const Export<
typename BlockCrsMatrixType::local_ordinal_type,
 
   37                                                 typename BlockCrsMatrixType::global_ordinal_type,
 
   38                                                 typename BlockCrsMatrixType::node_type>& exporter);
 
  120#if defined(TPETRA_ENABLE_BLOCKCRS_LITTLEBLOCK_LAYOUTLEFT) 
  136  using STS              = Teuchos::ScalarTraits<Scalar>;
 
  177  typedef ::Tpetra::MultiVector<Scalar, LO, GO, node_type> 
mv_type;
 
  185                       Kokkos::MemoryTraits<Kokkos::Unmanaged> >
 
  187  typedef typename little_block_type::host_mirror_type little_block_host_type;
 
  193                       Kokkos::MemoryTraits<Kokkos::Unmanaged> >
 
  197  typedef typename BMV::little_host_vec_type little_host_vec_type;
 
  201  typedef typename BMV::const_little_host_vec_type const_host_little_vec_type;
 
  204  using local_inds_device_view_type =
 
  205      typename row_matrix_type::local_inds_device_view_type;
 
  206  using local_inds_host_view_type =
 
  207      typename row_matrix_type::local_inds_host_view_type;
 
  208  using nonconst_local_inds_host_view_type =
 
  209      typename row_matrix_type::nonconst_local_inds_host_view_type;
 
  211  using global_inds_device_view_type =
 
  212      typename row_matrix_type::global_inds_device_view_type;
 
  213  using global_inds_host_view_type =
 
  214      typename row_matrix_type::global_inds_host_view_type;
 
  215  using nonconst_global_inds_host_view_type =
 
  216      typename row_matrix_type::nonconst_global_inds_host_view_type;
 
  218  using values_device_view_type =
 
  219      typename row_matrix_type::values_device_view_type;
 
  220  using values_host_view_type =
 
  221      typename row_matrix_type::values_host_view_type;
 
  222  using nonconst_values_host_view_type =
 
  223      typename row_matrix_type::nonconst_values_host_view_type;
 
  227  using local_matrix_device_type =
 
  232                                            typename local_graph_device_type::size_type>;
 
  237  using local_matrix_int_rowptrs_device_type =
 
  246      local_matrix_device_type,
 
  247      local_matrix_int_rowptrs_device_type,
 
  248      typename mv_type::device_view_type>;
 
  257  mutable std::shared_ptr<ApplyHelper> applyHelper;
 
  260  std::shared_ptr<ApplyHelper> getApplyHelper()
 const {
 
  263      applyHelper = std::make_shared<ApplyHelper>(A_lcl.nnz(), A_lcl.graph.row_map);
 
  268#if KOKKOS_VERSION >= 40799 
  269  using local_matrix_host_type =
 
  270      typename local_matrix_device_type::host_mirror_type;
 
  272  using local_matrix_host_type =
 
  273      typename local_matrix_device_type::HostMirror;
 
  295                 const typename local_matrix_device_type::values_type& 
values,
 
  318  Teuchos::RCP<const map_type> 
getDomainMap() 
const override;
 
  321  Teuchos::RCP<const map_type> 
getRangeMap() 
const override;
 
  324  Teuchos::RCP<const map_type> 
getRowMap() 
const override;
 
  327  Teuchos::RCP<const map_type> 
getColMap() 
const override;
 
  349        Teuchos::ETransp 
mode = Teuchos::NO_TRANS,
 
  350        Scalar alpha          = Teuchos::ScalarTraits<Scalar>::one(),
 
  351        Scalar beta           = Teuchos::ScalarTraits<Scalar>::zero()) 
const override;
 
  396           const Teuchos::EVerbosityLevel 
verbLevel) 
const override;
 
  406  virtual Teuchos::RCP<const ::Tpetra::RowGraph<LO, GO, Node> > 
getGraph() 
const override;
 
  415  applyBlock(
const BlockMultiVector<Scalar, LO, GO, Node>& X,
 
  416             BlockMultiVector<Scalar, LO, GO, Node>& Y,
 
  417             Teuchos::ETransp mode = Teuchos::NO_TRANS,
 
  418             const Scalar alpha    = Teuchos::ScalarTraits<Scalar>::one(),
 
  419             const Scalar beta     = Teuchos::ScalarTraits<Scalar>::zero());
 
  425                        const Import<LO, GO, Node>& importer) 
const;
 
  431                        const Export<LO, GO, Node>& exporter) 
const;
 
  462                        const LO numColInds) 
const;
 
  493                        const LO numColInds) 
const;
 
  529                  local_inds_host_view_type& indices,
 
  530                  values_host_view_type& values) 
const override;
 
  536                          local_inds_host_view_type& indices,
 
  537                          nonconst_values_host_view_type& values) 
const;
 
  542                  nonconst_local_inds_host_view_type& Indices,
 
  543                  nonconst_values_host_view_type& Values,
 
  544                  size_t& NumEntries) 
const override;
 
  546  getLocalBlockDeviceNonConst(
const LO localRowInd, 
const LO localColInd) 
const;
 
  548  little_block_host_type
 
  549  getLocalBlockHostNonConst(
const LO localRowInd, 
const LO localColInd) 
const;
 
  577                        const LO numColInds) 
const;
 
  585                                 const ptrdiff_t offsets[],
 
  587                                 const LO numOffsets) 
const;
 
  589  LO absMaxLocalValuesByOffsets(
const LO localRowInd,
 
  590                                const ptrdiff_t offsets[],
 
  592                                const LO numOffsets) 
const;
 
  600                                 const ptrdiff_t offsets[],
 
  602                                 const LO numOffsets) 
const;
 
  651    return (*errs_).is_null() ? std::string(
"") : (*errs_)->str();
 
 
  687                                         Kokkos::MemoryUnmanaged>& offsets) 
const;
 
  704                                      Kokkos::MemoryUnmanaged>& diag,
 
  706                                      Kokkos::MemoryUnmanaged>& offsets) 
const;
 
  742  virtual bool checkSizes(const ::Tpetra::SrcDistObject& 
source) 
override;
 
  750                 const size_t numSameIDs,
 
  768                 Kokkos::DualView<
size_t*,
 
  784                   Kokkos::DualView<
size_t*,
 
  794  Teuchos::RCP<crs_graph_type> graphRCP_;
 
  834  using graph_row_offset_host_type = 
typename crs_graph_type::local_graph_device_type::row_map_type::host_mirror_type;
 
  835  graph_row_offset_host_type ptrHost_;
 
  842  using graph_column_indices_host_type = 
typename crs_graph_type::local_graph_device_type::entries_type::host_mirror_type;
 
  843  graph_column_indices_host_type indHost_;
 
  850  using impl_scalar_type_dualview         = Kokkos::DualView<impl_scalar_type*, device_type>;
 
  875  Teuchos::RCP<Teuchos::RCP<BMV> > X_colMap_;
 
  879  Teuchos::RCP<Teuchos::RCP<BMV> > Y_rowMap_;
 
  888  Teuchos::RCP<Teuchos::RCP<typename crs_graph_type::import_type> > pointImporter_;
 
  904  Teuchos::RCP<bool> localError_;
 
  913  Teuchos::RCP<Teuchos::RCP<std::ostringstream> > errs_;
 
  916  std::ostream& markLocalErrorAndGetStream();
 
  922  typename impl_scalar_type_dualview::t_host::const_type
 
  923  getValuesHost() 
const;
 
  925  typename impl_scalar_type_dualview::t_dev::const_type
 
  926  getValuesDevice() 
const;
 
  946  typename impl_scalar_type_dualview::t_host
 
  949  typename impl_scalar_type_dualview::t_dev
 
  950  getValuesDeviceNonConst() 
const;
 
  953  typename impl_scalar_type_dualview::t_host::const_type
 
  954  getValuesHost(
const LO& 
lclRow) 
const;
 
  957  typename impl_scalar_type_dualview::t_dev::const_type
 
  958  getValuesDevice(
const LO& 
lclRow) 
const;
 
  961  typename impl_scalar_type_dualview::t_host
 
  965  typename impl_scalar_type_dualview::t_dev
 
  966  getValuesDeviceNonConst(
const LO& 
lclRow);
 
  982                  const Teuchos::ETransp 
mode,
 
 1053                                const LO 
hint = 0) 
const;
 
 1057  LO offsetPerBlock() 
const;
 
 1065  little_block_host_type
 
 1070  virtual Teuchos::RCP<const Teuchos::Comm<int> > 
getComm() 
const override;
 
 1100  virtual bool hasColMap() 
const override;
 
 1156                   nonconst_global_inds_host_view_type& 
Indices,
 
 1157                   nonconst_values_host_view_type& 
Values,
 
 1185                   global_inds_host_view_type& indices,
 
 1186                   values_host_view_type& 
values) 
const override;
 
 1210  virtual void leftScale(const ::Tpetra::Vector<Scalar, LO, GO, Node>& 
x) 
override;
 
 1217  virtual void rightScale(const ::Tpetra::Vector<Scalar, LO, GO, Node>& 
x) 
override;
 
 1227  virtual typename ::Tpetra::RowMatrix<Scalar, LO, GO, Node>::mag_type
 
 1232  template <
class BlockCrsMatrixType>
 
 1233  friend Teuchos::RCP<BlockCrsMatrixType>
 
 1234  Tpetra::importAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& 
sourceMatrix,
 
 1235                                              const Import<
typename BlockCrsMatrixType::local_ordinal_type,
 
 1236                                                           typename BlockCrsMatrixType::global_ordinal_type,
 
 1237                                                           typename BlockCrsMatrixType::node_type>& 
importer);
 
 1239  template <
class BlockCrsMatrixType>
 
 1240  friend Teuchos::RCP<BlockCrsMatrixType>
 
 1241  Tpetra::exportAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& 
sourceMatrix,
 
 1242                                              const Export<
typename BlockCrsMatrixType::local_ordinal_type,
 
 1243                                                           typename BlockCrsMatrixType::global_ordinal_type,
 
 1244                                                           typename BlockCrsMatrixType::node_type>& 
exporter);
 
 
 1247template <
class BlockCrsMatrixType>
 
 1248Teuchos::RCP<BlockCrsMatrixType>
 
 1249importAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& 
sourceMatrix,
 
 1250                                    const Import<
typename BlockCrsMatrixType::local_ordinal_type,
 
 1251                                                 typename BlockCrsMatrixType::global_ordinal_type,
 
 1252                                                 typename BlockCrsMatrixType::node_type>& 
importer) {
 
 1258template <
class BlockCrsMatrixType>
 
 1259Teuchos::RCP<BlockCrsMatrixType>
 
 1260exportAndFillCompleteBlockCrsMatrix(
const Teuchos::RCP<const BlockCrsMatrixType>& sourceMatrix,
 
 1261                                    const Export<
typename BlockCrsMatrixType::local_ordinal_type,
 
 1262                                                 typename BlockCrsMatrixType::global_ordinal_type,
 
 1263                                                 typename BlockCrsMatrixType::node_type>& exporter) {
 
 1264  Teuchos::RCP<BlockCrsMatrixType> destMatrix;
 
 1265  sourceMatrix->exportAndFillComplete(destMatrix, exporter);
 
Declaration of the Tpetra::CrsMatrix class.
 
Sparse matrix whose entries are small dense square blocks, all of the same dimensions.
 
LO sumIntoLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Sum into values at the given (mesh, i.e., block) column indices, in the given (mesh,...
 
virtual ~BlockCrsMatrix()
Destructor (declared virtual for memory safety).
 
virtual LO getBlockSize() const override
The number of degrees of freedom per mesh point.
 
void exportAndFillComplete(Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > &destMatrix, const Export< LO, GO, Node > &exporter) const
Import from this to the given destination matrix, and make the result fill complete.
 
virtual bool isLocallyIndexed() const override
Whether matrix indices are locally indexed.
 
virtual bool isFillComplete() const override
Whether fillComplete() has been called.
 
virtual size_t getGlobalMaxNumRowEntries() const override
The maximum number of entries in any row over all processes in the matrix's communicator.
 
void applyBlock(const BlockMultiVector< Scalar, LO, GO, Node > &X, BlockMultiVector< Scalar, LO, GO, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, const Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), const Scalar beta=Teuchos::ScalarTraits< Scalar >::zero())
Version of apply() that takes BlockMultiVector input and output.
 
LO replaceLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like replaceLocalValues, but avoids computing row offsets.
 
virtual void getLocalRowCopy(LO LocalRow, nonconst_local_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const override
Not implemented.
 
std::string errorMessages() const
The current stream of error messages.
 
LO absMaxLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos:...
 
virtual void getGlobalRowView(GO GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const override
Get a constant, nonpersisting, globally indexed view of the given row of the matrix.
 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
Print a description of this object to the given output stream.
 
size_t getLocalMaxNumRowEntries() const override
Maximum number of entries in any row of the matrix, on this process.
 
LO local_ordinal_type
The type of local indices.
 
virtual void getGlobalRowCopy(GO GlobalRow, nonconst_global_inds_host_view_type &Indices, nonconst_values_host_view_type &Values, size_t &NumEntries) const override
Get a copy of the given global row's entries.
 
LO getLocalRowOffsets(const LO localRowInd, ptrdiff_t offsets[], const LO colInds[], const LO numColInds) const
Get relative offsets corresponding to the given rows, given by local row index.
 
virtual size_t getLocalNumCols() const override
The number of columns needed to apply the forward operator on this node.
 
virtual void copyAndPermute(const SrcDistObject &sourceObj, const size_t numSameIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteToLIDs, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &permuteFromLIDs, const CombineMode CM) override
 
::Tpetra::MultiVector< Scalar, LO, GO, node_type > mv_type
The implementation of MultiVector that this class uses.
 
Teuchos::RCP< const map_type > getRangeMap() const override
Get the (point) range Map of this matrix.
 
LO sumIntoLocalValuesByOffsets(const LO localRowInd, const ptrdiff_t offsets[], const Scalar vals[], const LO numOffsets) const
Like sumIntoLocalValues, but avoids computing row offsets.
 
virtual bool hasColMap() const override
Whether this matrix has a well-defined column Map.
 
device_type::execution_space execution_space
The Kokkos execution space that this class uses.
 
size_t getLocalNumRows() const override
get the local number of block rows
 
virtual size_t getLocalNumEntries() const override
The local number of stored (structurally nonzero) entries.
 
bool hasTransposeApply() const override
Whether it is valid to apply the transpose or conjugate transpose of this matrix.
 
impl_scalar_type_dualview::t_host getValuesHostNonConst() const
Get the host or device View of the matrix's values (val_).
 
Kokkos::View< impl_scalar_type **, Impl::BlockCrsMatrixLittleBlockArrayLayout, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > little_block_type
The type used to access nonconst matrix blocks.
 
Teuchos::RCP< const map_type > getColMap() const override
get the (mesh) map for the columns of this block matrix.
 
virtual typename::Tpetra::RowMatrix< Scalar, LO, GO, Node >::mag_type getFrobeniusNorm() const override
The Frobenius norm of the matrix.
 
virtual Teuchos::RCP< const ::Tpetra::RowGraph< LO, GO, Node > > getGraph() const override
Get the (mesh) graph.
 
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode) override
 
local_matrix_device_type getLocalMatrixDevice() const
 
void apply(const mv_type &X, mv_type &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=Teuchos::ScalarTraits< Scalar >::one(), Scalar beta=Teuchos::ScalarTraits< Scalar >::zero()) const override
For this matrix A, compute Y := beta * Y + alpha * Op(A) * X.
 
virtual global_size_t getGlobalNumCols() const override
The global number of columns of this matrix.
 
void getLocalDiagOffsets(const Kokkos::View< size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Get offsets of the diagonal entries in the matrix.
 
Teuchos::RCP< const map_type > getRowMap() const override
get the (mesh) map for the rows of this block matrix.
 
std::string description() const override
One-line description of this object.
 
void importAndFillComplete(Teuchos::RCP< BlockCrsMatrix< Scalar, LO, GO, Node > > &destMatrix, const Import< LO, GO, Node > &importer) const
Import from this to the given destination matrix, and make the result fill complete.
 
Node node_type
The Node type.
 
void getLocalRowView(LO LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const override
Get a view of the (mesh, i.e., block) row, using local (mesh, i.e., block) indices.
 
void getLocalDiagCopy(const Kokkos::View< impl_scalar_type ***, device_type, Kokkos::MemoryUnmanaged > &diag, const Kokkos::View< const size_t *, device_type, Kokkos::MemoryUnmanaged > &offsets) const
Variant of getLocalDiagCopy() that uses precomputed offsets and puts diagonal blocks in a 3-D Kokkos:...
 
GO global_ordinal_type
The type of global indices.
 
size_t getNumEntriesInLocalRow(const LO localRowInd) const override
Return the number of entries in the given row on the calling process.
 
virtual bool supportsRowViews() const override
Whether this object implements getLocalRowView() and getGlobalRowView().
 
LO replaceLocalValues(const LO localRowInd, const LO colInds[], const Scalar vals[], const LO numColInds) const
Replace values at the given (mesh, i.e., block) column indices, in the given (mesh,...
 
virtual size_t getNumEntriesInGlobalRow(GO globalRow) const override
The current number of entries on the calling process in the specified global row.
 
Node::device_type device_type
The Kokkos::Device specialization that this class uses.
 
BMV::little_vec_type little_vec_type
The type used to access nonconst vector blocks.
 
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to alpha.
 
BMV::const_little_vec_type const_little_vec_type
The type used to access const vector blocks.
 
Kokkos::View< const impl_scalar_type **, Impl::BlockCrsMatrixLittleBlockArrayLayout, device_type, Kokkos::MemoryTraits< Kokkos::Unmanaged > > const_little_block_type
The type used to access const matrix blocks.
 
void getLocalRowViewNonConst(LO LocalRow, local_inds_host_view_type &indices, nonconst_values_host_view_type &values) const
 
char packet_type
Implementation detail; tells.
 
typename DistObject< Scalar, LO, GO, Node >::buffer_device_type buffer_device_type
Kokkos::Device specialization for communication buffers.
 
virtual void packAndPrepare(const SrcDistObject &sourceObj, 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) override
 
::Tpetra::Map< LO, GO, node_type > map_type
The implementation of Map that this class uses.
 
virtual void rightScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x) override
Scale the RowMatrix on the right with the given Vector x.
 
virtual GO getIndexBase() const override
The index base for global indices in this matrix.
 
bool localError() const
Whether this object had an error on the calling process.
 
virtual global_size_t getGlobalNumEntries() const override
The global number of stored (structurally nonzero) entries.
 
typename BMV::impl_scalar_type impl_scalar_type
The implementation type of entries in the matrix.
 
virtual Teuchos::RCP< const Teuchos::Comm< int > > getComm() const override
The communicator over which this matrix is distributed.
 
::Tpetra::CrsGraph< LO, GO, node_type > crs_graph_type
The implementation of CrsGraph that this class uses.
 
device_type::memory_space memory_space
The Kokkos memory space that this class uses.
 
global_size_t getGlobalNumRows() const override
get the global number of block rows
 
virtual void leftScale(const ::Tpetra::Vector< Scalar, LO, GO, Node > &x) override
Scale the RowMatrix on the left with the given Vector x.
 
Teuchos::RCP< const map_type > getDomainMap() const override
Get the (point) domain Map of this matrix.
 
virtual bool isGloballyIndexed() const override
Whether matrix indices are globally indexed.
 
BlockCrsMatrix()
Default constructor: Makes an empty block matrix.
 
MultiVector for multiple degrees of freedom per mesh point.
 
Kokkos::View< impl_scalar_type *, device_type > little_vec_type
"Block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
 
typename mv_type::impl_scalar_type impl_scalar_type
The implementation type of entries in the object.
 
Kokkos::View< const impl_scalar_type *, device_type > const_little_vec_type
"Const block view" of all degrees of freedom at a mesh point, for a single column of the MultiVector.
 
KokkosSparse::StaticCrsGraph< local_ordinal_type, Kokkos::LayoutLeft, device_type, void, size_t > local_graph_device_type
The type of the part of the sparse graph on each MPI process.
 
Struct that holds views of the contents of a CrsMatrix.
 
Base class for distributed Tpetra objects that support data redistribution.
 
virtual void packAndPrepare(const SrcDistObject &source, const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &exportLIDs, Kokkos::DualView< packet_type *, buffer_device_type > &exports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, size_t &constantNumPackets)
Pack data and metadata for communication (sends).
 
virtual void unpackAndCombine(const Kokkos::DualView< const local_ordinal_type *, buffer_device_type > &importLIDs, Kokkos::DualView< packet_type *, buffer_device_type > imports, Kokkos::DualView< size_t *, buffer_device_type > numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode)
Perform any unpacking and combining after communication.
 
Communication plan for data redistribution from a (possibly) multiply-owned to a uniquely-owned distr...
 
Communication plan for data redistribution from a uniquely-owned to a (possibly) multiply-owned distr...
 
A read-only, row-oriented interface to a sparse matrix.
 
Abstract base class for objects that can be the source of an Import or Export operation.
 
Kokkos::LayoutRight BlockCrsMatrixLittleBlockArrayLayout
give an option to use layoutleft
 
Namespace Tpetra contains the class and methods constituting the Tpetra library.
 
size_t global_size_t
Global size_t object.
 
CombineMode
Rule for combining data in an Import or Export.