10 #ifndef XPETRA_TPETRACRSMATRIX_DECL_HPP 
   11 #define XPETRA_TPETRACRSMATRIX_DECL_HPP 
   21 #include "Tpetra_CrsMatrix.hpp" 
   22 #include "Tpetra_FECrsMatrix.hpp" 
   23 #include "Tpetra_replaceDiagonalCrsMatrix.hpp" 
   31 #include "Xpetra_TpetraCrsGraph.hpp" 
   36 template <
class Scalar,
 
   39           class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
 
   41   : 
public CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>  
 
   43 #undef XPETRA_TPETRACRSMATRIX_SHORT 
   53 #ifdef HAVE_XPETRA_TPETRA 
   84                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null);
 
   91                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null);
 
   99                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms);
 
  107                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms);
 
  109 #ifdef HAVE_XPETRA_TPETRA 
  133                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null);
 
  142       const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = null);
 
  153       const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null);
 
  165   void insertGlobalValues(GlobalOrdinal globalRow, 
const ArrayView<const GlobalOrdinal> &cols, 
const ArrayView<const Scalar> &vals);
 
  168   void insertLocalValues(LocalOrdinal localRow, 
const ArrayView<const LocalOrdinal> &cols, 
const ArrayView<const Scalar> &vals);
 
  171   void replaceGlobalValues(GlobalOrdinal globalRow, 
const ArrayView<const GlobalOrdinal> &cols, 
const ArrayView<const Scalar> &vals);
 
  176                      const ArrayView<const LocalOrdinal> &cols,
 
  177                      const ArrayView<const Scalar> &vals);
 
  183   void scale(
const Scalar &alpha);
 
  187   void allocateAllValues(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values);
 
  190   void setAllValues(
const ArrayRCP<size_t> &rowptr, 
const ArrayRCP<LocalOrdinal> &colind, 
const ArrayRCP<Scalar> &values);
 
  193   void getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values) 
const;
 
  206   void resumeFill(
const RCP<ParameterList> ¶ms = null);
 
  212   void fillComplete(
const RCP<ParameterList> ¶ms = null);
 
  222                                 const RCP<ParameterList> ¶ms                                     = Teuchos::null);
 
  230   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRowMap() 
const;
 
  233   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getColMap() 
const;
 
  236   RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> 
getCrsGraph() 
const;
 
  287   void getLocalRowCopy(LocalOrdinal LocalRow, 
const ArrayView<LocalOrdinal> &Indices, 
const ArrayView<Scalar> &Values, 
size_t &NumEntries) 
const;
 
  290   void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) 
const;
 
  293   void getGlobalRowCopy(GlobalOrdinal GlobalRow, 
const ArrayView<GlobalOrdinal> &indices, 
const ArrayView<Scalar> &values, 
size_t &numEntries) 
const;
 
  296   void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) 
const;
 
  304   void apply(
const MultiVector &X, 
MultiVector &Y, Teuchos::ETransp mode = Teuchos::NO_TRANS, Scalar alpha = ScalarTraits<Scalar>::one(), Scalar beta = ScalarTraits<Scalar>::zero()) 
const;
 
  310   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getDomainMap() 
const;
 
  313   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRangeMap() 
const;
 
  324   void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) 
const;
 
  346   void getLocalDiagCopy(
Vector &diag, 
const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) 
const;
 
  361   Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getMap() 
const;
 
  389   TpetraCrsMatrix(
const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx);
 
  392   RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> 
getTpetra_CrsMatrix() 
const;
 
  397 #ifdef HAVE_XPETRA_TPETRA 
  398 #if KOKKOS_VERSION >= 40799 
  412                     const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind,
 
  413                     const typename local_matrix_type::values_type &val) {
 
  429   RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> 
mtx_;
 
  432 #ifdef HAVE_XPETRA_EPETRA 
  434 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_INT))) || \ 
  435      (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_INT)))) 
  438 template <
class Scalar>
 
  440   : 
public CrsMatrix<Scalar, int, int, EpetraNode>  
 
  452 #ifdef HAVE_XPETRA_TPETRA 
  490                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null) {
 
  499                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null) {
 
  509                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
 
  519                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
 
  523 #ifdef HAVE_XPETRA_TPETRA 
  547                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null) {
 
  558       const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = null) {
 
  583                      const ArrayView<const LocalOrdinal> &cols,
 
  584                      const ArrayView<const Scalar> &vals) {}
 
  594   void allocateAllValues(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) {}
 
  597   void setAllValues(
const ArrayRCP<size_t> &rowptr, 
const ArrayRCP<LocalOrdinal> &colind, 
const ArrayRCP<Scalar> &values) {}
 
  600   void getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values)
 const {}
 
  629                                 const RCP<ParameterList> ¶ms                                     = Teuchos::null) {}
 
  637   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRowMap()
 const { 
return Teuchos::null; }
 
  640   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getColMap()
 const { 
return Teuchos::null; }
 
  643   RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> 
getCrsGraph()
 const { 
return Teuchos::null; }
 
  688   typename ScalarTraits<Scalar>::magnitudeType 
getFrobeniusNorm()
 const { 
return ScalarTraits<Scalar>::magnitude(ScalarTraits<Scalar>::zero()); }
 
  711   void apply(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, 
MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode = Teuchos::NO_TRANS, Scalar alpha = ScalarTraits<Scalar>::one(), Scalar beta = ScalarTraits<Scalar>::zero())
 const {}
 
  714   void apply(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, 
MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, 
bool sumInterfaceValues, 
const RCP<
Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> ®ionInterfaceImporter, 
const Teuchos::ArrayRCP<LocalOrdinal> ®ionInterfaceLIDs)
 const {}
 
  717   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getDomainMap()
 const { 
return Teuchos::null; }
 
  720   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRangeMap()
 const { 
return Teuchos::null; }
 
  731   void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
 const {}
 
  765   Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getMap()
 const { 
return Teuchos::null; }
 
  793   TpetraCrsMatrix(
const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx) {
 
  798   RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> 
getTpetra_CrsMatrix()
 const { 
return Teuchos::null; }
 
  803 #ifdef HAVE_XPETRA_TPETRA 
  810                     const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind,
 
  811                     const typename local_matrix_type::values_type &val) {}
 
  826 #if ((defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_OPENMP) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG))) || \ 
  827      (!defined(EPETRA_HAVE_OMP) && (!defined(HAVE_TPETRA_INST_SERIAL) || !defined(HAVE_TPETRA_INST_INT_LONG_LONG)))) 
  830 template <
class Scalar>
 
  832   : 
public CrsMatrix<Scalar, int, long long, EpetraNode>  
 
  844 #ifdef HAVE_XPETRA_TPETRA 
  882                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null) {
 
  891                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = Teuchos::null) {
 
  901                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
 
  911                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms) {
 
  915 #ifdef HAVE_XPETRA_TPETRA 
  939                   const Teuchos::RCP<Teuchos::ParameterList> ¶ms = null) {
 
  950       const Teuchos::RCP<Teuchos::ParameterList> ¶ms                          = null) {
 
  975                      const ArrayView<const LocalOrdinal> &cols,
 
  976                      const ArrayView<const Scalar> &vals) {}
 
  986   void allocateAllValues(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) {}
 
  989   void setAllValues(
const ArrayRCP<size_t> &rowptr, 
const ArrayRCP<LocalOrdinal> &colind, 
const ArrayRCP<Scalar> &values) {}
 
  992   void getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values)
 const {}
 
 1021                                 const RCP<ParameterList> ¶ms                                     = Teuchos::null) {}
 
 1029   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRowMap()
 const { 
return Teuchos::null; }
 
 1032   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getColMap()
 const { 
return Teuchos::null; }
 
 1035   RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node>> 
getCrsGraph()
 const { 
return Teuchos::null; }
 
 1080   typename ScalarTraits<Scalar>::magnitudeType 
getFrobeniusNorm()
 const { 
return ScalarTraits<Scalar>::magnitude(ScalarTraits<Scalar>::zero()); }
 
 1103   void apply(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, 
MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode = Teuchos::NO_TRANS, Scalar alpha = ScalarTraits<Scalar>::one(), Scalar beta = ScalarTraits<Scalar>::zero())
 const {}
 
 1106   void apply(
const MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &X, 
MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, 
bool sumInterfaceValues, 
const RCP<
Xpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> ®ionInterfaceImporter, 
const Teuchos::ArrayRCP<LocalOrdinal> ®ionInterfaceLIDs)
 const {}
 
 1109   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getDomainMap()
 const { 
return Teuchos::null; }
 
 1112   const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getRangeMap()
 const { 
return Teuchos::null; }
 
 1123   void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
 const {}
 
 1157   Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> 
getMap()
 const { 
return Teuchos::null; }
 
 1185   TpetraCrsMatrix(
const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx) {
 
 1190   RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> 
getTpetra_CrsMatrix()
 const { 
return Teuchos::null; }
 
 1195 #ifdef HAVE_XPETRA_TPETRA 
 1202 #if KOKKOS_VERSION >= 40799 
 1204     TEUCHOS_UNREACHABLE_RETURN(
typename local_matrix_type::host_mirror_type());
 
 1208     TEUCHOS_UNREACHABLE_RETURN(
typename local_matrix_type::HostMirror());
 
 1217                     const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind,
 
 1218                     const typename local_matrix_type::values_type &val) {}
 
 1233 #endif  // HAVE_XPETRA_EPETRA 
 1235 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1236 Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
 
 1241 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1242 Teuchos::RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>>
 
 1247 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1248 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
 
 1253 template <
class Scalar, 
class LocalOrdinal, 
class GlobalOrdinal, 
class Node>
 
 1254 Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> &
 
 1261 #define XPETRA_TPETRACRSMATRIX_SHORT 
 1262 #endif  // XPETRA_TPETRACRSMATRIX_DECL_HPP 
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import. 
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const 
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
LocalOrdinal GetStorageBlockSize() const 
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
size_t getGlobalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on all nodes. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const 
Returns the Map that describes the row distribution in this matrix. 
global_size_t getGlobalNumRows() const 
Number of global elements in the row map of this matrix. 
void replaceDiag(const Vector &diag)
Replace the diagonal entries of the matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row. 
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const 
Get the underlying Tpetra matrix. 
std::string description() const 
A simple one-line description of this object. 
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const 
Returns the current number of entries on this node in the specified local row. 
void fillComplete(const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete. 
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import. 
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis. 
void leftScale(const Vector &x)
Left scale operator with given vector values. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const 
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const 
Implements DistObject interface. 
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs. 
global_size_t getGlobalNumCols() const 
Number of global columns in the matrix. 
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &newImporter)
Replaces the current domainMap and importer with the user-specified objects. 
LocalOrdinal GetStorageBlockSize() const 
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export. 
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const 
Returns the Map that describes the column distribution in this matrix. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const 
Returns the Map associated with the range of this operator, which must be compatible with Y...
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row. 
bool isLocallyIndexed() const 
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
global_size_t getGlobalNumEntries() const 
Returns the global number of entries in this matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node >> DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor for a fused export. 
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mtx_
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node >> &exporter=Teuchos::null, const RCP< ParameterList > ¶ms=Teuchos::null)
Expert static fill complete. 
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps. 
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const 
Get offsets of the diagonal entries in the matrix. 
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of local indices in a specified row of the matrix. 
void getAllValues(ArrayRCP< Scalar > &values)
Gets the 1D pointer arrays of the graph. 
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const 
Returns the current number of entries in the (locally owned) global row. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying fixed number of entries for each row. 
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const 
Get a copy of the diagonal entries owned by this node, with local row idices. 
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const 
Compute a residual R = B - (*this) * X. 
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph. 
void setAllValues(const typename local_matrix_type::row_map_type &ptr, const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind, const typename local_matrix_type::values_type &val)
size_t getLocalNumRows() const 
Returns the number of matrix rows owned on the calling node. 
size_t getLocalNumEntries() const 
Returns the local number of entries in this matrix. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const 
Returns the Map associated with the range of this operator, which must be compatible with Y...
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values. 
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const 
Returns the Frobenius norm of the matrix. 
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine. 
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs. 
void resumeFill(const RCP< ParameterList > ¶ms=null)
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs. 
bool isLocallyIndexed() const 
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void setObjectLabel(const std::string &objectLabel)
void setObjectLabel(const std::string &objectLabel)
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of local indices in a specified row of the matrix. 
TpetraVector TpetraVectorClass
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
void setObjectLabel(const std::string &objectLabel)
local_matrix_type::HostMirror getLocalMatrixHost() const 
Access the local Kokkos::CrsMatrix data. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &RowExporter, const Teuchos::RCP< const Export< LocalOrdinal, GlobalOrdinal, Node >> DomainExporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor for a fused export. 
bool haveGlobalConstants() const 
Returns true if globalConstants have been computed; false otherwise. 
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newMap)
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph. 
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const 
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
virtual ~TpetraCrsMatrix()
Destructor. 
global_size_t getGlobalNumRows() const 
Number of global elements in the row map of this matrix. 
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const 
Returns the CrsGraph associated with this matrix. 
void getAllValues(ArrayRCP< Scalar > &values)
Gets the 1D pointer arrays of the graph. 
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer). 
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this. 
#define XPETRA_TPETRA_ETI_EXCEPTION(cl, obj, go, node)
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newMap)
local_matrix_type::HostMirror getLocalMatrixHost() const 
Access the local Kokkos::CrsMatrix data. 
local_matrix_type getLocalMatrixDevice() const 
Access the local Kokkos::CrsMatrix data. 
void fillComplete(const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete. 
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const 
Get a copy of the diagonal entries owned by this node, with local row indices. 
void getLocalDiagCopy(Vector &diag) const 
Get a copy of the diagonal entries owned by this node, with local row indices. 
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps. 
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values. 
std::string description() const 
A simple one-line description of this object. 
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const 
Gets the 1D pointer arrays of the graph. 
size_t getLocalNumCols() const 
Returns the number of columns connected to the locally owned rows of this matrix. ...
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor for a fused export. 
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const 
Returns the Map that describes the row distribution in this matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node >> &graph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying a previously constructed graph. 
bool isFillComplete() const 
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called. 
bool isLocallyIndexed() const 
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis. 
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node >> &graph, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying a previously constructed graph. 
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const 
Get offsets of the diagonal entries in the matrix. 
std::string description() const 
A simple one-line description of this object. 
bool isFillComplete() const 
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying fixed number of entries for each row. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const 
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node >> ®ionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > ®ionInterfaceLIDs) const 
Computes the matrix-multivector multiplication for region layout matrices. 
TpetraCrsMatrix(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &mtx)
TpetraCrsMatrix constructor to wrap a Tpetra::CrsMatrix object. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column Map and number of entries in each row. 
bool haveGlobalConstants() const 
Returns true if globalConstants have been computed; false otherwise. 
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const 
Compute a residual R = B - (*this) * X. 
size_t getLocalNumRows() const 
Returns the number of matrix rows owned on the calling node. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor for a fused import. 
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const 
Get the underlying Tpetra matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying fixed number of entries for each row. 
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const 
Implements DistObject interface. 
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const 
Get offsets of the diagonal entries in the matrix. 
bool isGloballyIndexed() const 
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void setAllValues(const typename local_matrix_type::row_map_type &ptr, const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind, const typename local_matrix_type::values_type &val)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const 
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node >> &exporter=Teuchos::null, const RCP< ParameterList > ¶ms=Teuchos::null)
Expert static fill complete. 
global_size_t getGlobalNumRows() const 
Number of global elements in the row map of this matrix. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const 
Returns the Map that describes the column distribution in this matrix. 
TpetraCrsMatrix(const local_matrix_type &lclMatrix, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Constructor specifying local matrix and 4 maps. 
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const 
Implements DistObject interface. 
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const 
Returns the current number of entries in the (locally owned) global row. 
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const 
Returns the CrsGraph associated with this matrix. 
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const 
Compute a residual R = B - (*this) * X. 
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export. 
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const 
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter). 
TpetraCrsMatrix(const local_matrix_type &lclMatrix, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=null)
Constructor specifying local matrix and 4 maps. 
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const 
Computes the sparse matrix-multivector multiplication. 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Print the object with some verbosity level to an FancyOStream object. 
size_t getLocalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on this node. 
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine. 
size_t getLocalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on this node. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column Map and fixed number of entries for each row. 
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of global indices in a specified row of the matrix. 
void resumeFill(const RCP< ParameterList > ¶ms=null)
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Print the object with some verbosity level to an FancyOStream object. 
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this. 
size_t getGlobalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on all nodes. 
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of global indices in a specified row of the matrix. 
TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrixClass
size_t global_size_t
Global size_t object. 
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer). 
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs. 
global_size_t getGlobalNumCols() const 
Number of global columns in the matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor for a fused export. 
virtual ~TpetraCrsMatrix()
Destructor. 
TpetraVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraVectorClass
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const 
Returns the Map associated with the range of this operator, which must be compatible with Y...
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &RowImporter, const Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> DomainImporter, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Constructor for a fused import. 
bool hasMatrix() const 
Does this have an underlying matrix. 
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor for a fused import. 
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps. 
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const 
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
bool isGloballyIndexed() const 
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of global indices in a specified row of the matrix. 
virtual ~TpetraCrsMatrix()
Destructor. 
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph. 
size_t getLocalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on this node. 
bool isFillComplete() const 
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called. 
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Print the object with some verbosity level to an FancyOStream object. 
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs. 
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &newImporter)
Replaces the current domainMap and importer with the user-specified objects. 
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values. 
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const 
Get the underlying Tpetra matrix. 
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const 
Gets the 1D pointer arrays of the graph. 
void resumeFill(const RCP< ParameterList > ¶ms=null)
bool hasMatrix() const 
Does this have an underlying matrix. 
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const 
Computes the sparse matrix-multivector multiplication. 
global_size_t getGlobalNumEntries() const 
Returns the global number of entries in this matrix. 
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const 
Returns the current number of entries in the (locally owned) global row. 
bool haveGlobalConstants() const 
Returns true if globalConstants have been computed; false otherwise. 
global_size_t getGlobalNumCols() const 
Number of global columns in the matrix. 
bool supportsRowViews() const 
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class. 
global_size_t getGlobalNumEntries() const 
Returns the global number of entries in this matrix. 
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this. 
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const 
Get the underlying Tpetra matrix. 
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying column Map and number of entries in each row. 
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node >> ®ionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > ®ionInterfaceLIDs) const 
Computes the matrix-multivector multiplication for region layout matrices. 
size_t getLocalNumRows() const 
Returns the number of matrix rows owned on the calling node. 
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix. 
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &newImporter)
Replaces the current domainMap and importer with the user-specified objects. 
bool hasMatrix() const 
Does this have an underlying matrix. 
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const 
Returns the Map that describes the row distribution in this matrix. 
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
bool supportsRowViews() const 
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class. 
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying (possibly different) number of entries in each row. 
CombineMode
Xpetra::Combine Mode enumerable type. 
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor. 
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const 
Computes the sparse matrix-multivector multiplication. 
size_t getLocalNumCols() const 
Returns the number of columns connected to the locally owned rows of this matrix. ...
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const 
Get the underlying Tpetra matrix. 
size_t getLocalNumEntries() const 
Returns the local number of entries in this matrix. 
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const 
Get a copy of the diagonal entries owned by this node, with local row idices. 
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const 
Extract a const, non-persisting view of local indices in a specified row of the matrix. 
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export. 
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const 
Gets the 1D pointer arrays of the graph. 
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs. 
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine. 
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
void setAllValues(const typename local_matrix_type::row_map_type &ptr, const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind, const typename local_matrix_type::values_type &val)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const 
Returns the Map that describes the column distribution in this matrix. 
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values. 
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import. 
TpetraExport< LocalOrdinal, GlobalOrdinal, Node > TpetraExportClass
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &newMap)
bool isGloballyIndexed() const 
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix. 
bool isFillActive() const 
Returns true if the matrix is in edit mode. 
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const 
Returns the CrsGraph associated with this matrix. 
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const 
Returns the current number of entries on this node in the specified local row. 
size_t getGlobalMaxNumRowEntries() const 
Returns the maximum number of entries across all rows/columns on all nodes. 
LocalOrdinal GetStorageBlockSize() const 
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
size_t getLocalNumCols() const 
Returns the number of columns connected to the locally owned rows of this matrix. ...
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const 
Returns the Frobenius norm of the matrix. 
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const 
Returns the Frobenius norm of the matrix. 
bool isFillActive() const 
Returns true if the matrix is in edit mode. 
TpetraCrsMatrix(const TpetraCrsMatrix &matrix)
Deep copy constructor. 
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const 
Get the underlying Tpetra matrix. 
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs. 
TpetraImport< LocalOrdinal, GlobalOrdinal, Node > TpetraImportClass
TpetraCrsMatrix(const Teuchos::RCP< const CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &sourceMatrix, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap=Teuchos::null, const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap=Teuchos::null, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor for a fused import. 
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const 
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const 
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis. 
size_t getLocalNumEntries() const 
Returns the local number of entries in this matrix. 
bool isFillActive() const 
Returns true if the matrix is in edit mode. 
bool supportsRowViews() const 
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class. 
void expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node >> &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node >> &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node >> &exporter=Teuchos::null, const RCP< ParameterList > ¶ms=Teuchos::null)
Expert static fill complete. 
void rightScale(const Vector &x)
Right scale operator with given vector values. 
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const 
Returns the current number of entries on this node in the specified local row. 
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter). 
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const 
Get a copy of the diagonal entries owned by this node, with local row indices.