46 #ifndef XPETRA_CRSMATRIX_HPP 
   47 #define XPETRA_CRSMATRIX_HPP 
   58 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
   59 #ifdef HAVE_XPETRA_TPETRA 
   60 #include <Kokkos_StaticCrsGraph.hpp> 
   61 #include <KokkosSparse_CrsMatrix.hpp> 
   67   template <
class Scalar,
 
   72     : 
public RowMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>,
 
   73       public DistObject<char, LocalOrdinal,GlobalOrdinal,Node>
 
   94                         const ArrayView<const GlobalOrdinal>& cols,
 
   95                         const ArrayView<const Scalar>& vals) = 0;
 
   98     virtual void insertLocalValues(LocalOrdinal localRow, 
const ArrayView< const LocalOrdinal > &cols, 
const ArrayView< const Scalar > &vals)= 0;
 
  101     virtual void replaceGlobalValues(GlobalOrdinal globalRow, 
const ArrayView< const GlobalOrdinal > &cols, 
const ArrayView< const Scalar > &vals)= 0;
 
  104     virtual void replaceLocalValues(LocalOrdinal localRow, 
const ArrayView< const LocalOrdinal > &cols, 
const ArrayView< const Scalar > &vals)= 0;
 
  110     virtual void scale(
const Scalar &alpha)= 0;
 
  122     virtual void allocateAllValues(
size_t numNonZeros,ArrayRCP<size_t> & rowptr, ArrayRCP<LocalOrdinal> & colind, ArrayRCP<Scalar> & values)=0;
 
  125     virtual void setAllValues(
const ArrayRCP<size_t> & rowptr, 
const ArrayRCP<LocalOrdinal> & colind, 
const ArrayRCP<Scalar> & values)=0;
 
  128     virtual void getAllValues(ArrayRCP<const size_t>& rowptr, ArrayRCP<const LocalOrdinal>& colind, ArrayRCP<const Scalar>& values) 
const = 0;
 
  136     virtual void resumeFill(
const RCP< ParameterList > ¶ms=null)= 0;
 
  142     virtual void fillComplete(
const RCP< ParameterList > ¶ms=null)= 0;
 
  152                                           const RCP<ParameterList> ¶ms=Teuchos::null) = 0;
 
  159     virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRowMap() 
const = 0;
 
  162     virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getColMap() 
const = 0;
 
  165     virtual RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > 
getCrsGraph() 
const = 0;
 
  207     virtual typename ScalarTraits< Scalar >::magnitudeType 
getFrobeniusNorm() 
const = 0;
 
  213     virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) 
const = 0;
 
  216     virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, 
const ArrayView<GlobalOrdinal> &indices, 
const ArrayView<Scalar> &values, 
size_t &numEntries) 
const = 0;
 
  219     virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) 
const = 0;
 
  261     virtual 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 = 0;
 
  264     virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getDomainMap() 
const = 0;
 
  267     virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > >  
getRangeMap() 
const = 0;
 
  278     virtual void describe(Teuchos::FancyOStream &out, 
const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) 
const = 0;
 
  289 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR 
  290 #ifdef HAVE_XPETRA_TPETRA 
  291     typedef typename Kokkos::Details::ArithTraits<Scalar>::val_type impl_scalar_type;
 
  292     typedef typename node_type::execution_space execution_space;
 
  295     typedef Kokkos::StaticCrsGraph<LocalOrdinal,
 
  297                                        execution_space> local_graph_type;
 
  301     typedef KokkosSparse::CrsMatrix<impl_scalar_type, LocalOrdinal, execution_space,void,
 
  302                               typename local_graph_type::size_type> local_matrix_type;
 
  305     virtual local_matrix_type getLocalMatrix () 
const = 0;
 
  307     virtual void setAllValues (
const typename local_matrix_type::row_map_type& ptr,
 
  308                                const typename local_graph_type::entries_type::non_const_type& ind,
 
  309                                const typename local_matrix_type::values_type& val)=0;
 
  312 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too." 
  325         virtual void getLocalRowCopy(LocalOrdinal LocalRow, 
const ArrayView< LocalOrdinal > &Indices, 
const ArrayView< Scalar > &Values, 
size_t &NumEntries) 
const = 0;
 
  340 #define XPETRA_CRSMATRIX_SHORT 
  341 #endif // XPETRA_CRSMATRIX_HPP 
virtual 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)=0
Expert static fill complete. 
virtual global_size_t getGlobalNumRows() const =0
Number of global elements in the row map of this matrix. 
virtual 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 =0
Computes the sparse matrix-multivector multiplication. 
virtual size_t getNodeMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node. 
virtual void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)=0
Replace the diagonal entries of the matrix. 
virtual void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const =0
Compute a residual R = B - (*this) * X. 
virtual void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using local IDs. 
virtual void removeEmptyProcessesInPlace(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)=0
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
Returns the Map that describes the column distribution in this matrix. 
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Right scale matrix using the given vector entries. 
GlobalOrdinal global_ordinal_type
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this. 
virtual bool hasMatrix() const =0
Does this have an underlying matrix. 
virtual void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null)=0
Signal that data entry is complete, specifying domain and range maps. 
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs. 
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
Returns the current number of entries in the specified global row. 
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const =0
Returns the Map associated with the range of this operator, which must be compatible with Y...
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
Returns the Map that describes the row distribution in this matrix. 
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
Returns the current number of entries on this node in the specified local row. 
virtual bool isLocallyIndexed() const =0
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
virtual void setAllToScalar(const Scalar &alpha)=0
Set all matrix entries equal to scalarThis. 
LocalOrdinal local_ordinal_type
virtual void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)=0
Sets the 1D pointer arrays of the graph. 
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const =0
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
virtual size_t getNodeNumEntries() const =0
Returns the local number of entries in this matrix. 
virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const =0
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
virtual std::string description() const =0
A simple one-line description of this object. 
virtual void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)=0
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine. 
virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const =0
Extract a const, non-persisting view of local indices in a specified row of the matrix. 
virtual global_size_t getGlobalNumCols() const =0
Number of global columns in the matrix. 
virtual bool supportsRowViews() const =0
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class. 
virtual void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using global IDs. 
virtual void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const =0
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const =0
Print the object with some verbosity level to an FancyOStream object. 
virtual void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)=0
Replaces the current domainMap and importer with the user-specified objects. 
virtual void resumeFill(const RCP< ParameterList > ¶ms=null)=0
virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using local IDs. 
virtual bool haveGlobalConstants() const =0
Returns true if globalConstants have been computed; false otherwise. 
size_t global_size_t
Global size_t object. 
virtual global_size_t getGlobalNumEntries() const =0
Returns the global number of entries in this matrix. 
virtual size_t getGlobalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on all nodes. 
virtual void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const =0
Get a copy of the diagonal entries owned by this node, with local row indices. 
virtual ~CrsMatrix()
Destructor. 
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const =0
Returns the Frobenius norm of the matrix. 
virtual void setObjectLabel(const std::string &objectLabel)=0
virtual bool isFillComplete() const =0
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called. 
virtual void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const =0
Get offsets of the diagonal entries in the matrix. 
virtual size_t getNodeNumCols() const =0
Returns the number of matrix columns owned on the calling node. 
virtual RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const =0
Returns the CrsGraph associated with this matrix. 
virtual bool isGloballyIndexed() const =0
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
virtual size_t getNodeNumRows() const =0
Returns the number of matrix rows owned on the calling node. 
virtual bool isFillActive() const =0
Returns true if the matrix is in edit mode. 
virtual void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const =0
Gets the 1D pointer arrays of the graph. 
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Left scale matrix using the given vector entries. 
virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const =0
Extract a const, non-persisting view of global indices in a specified row of the matrix.