10 #ifndef XPETRA_CRSMATRIX_HPP
11 #define XPETRA_CRSMATRIX_HPP
15 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
22 #ifdef HAVE_XPETRA_TPETRA
23 #include <Kokkos_StaticCrsGraph.hpp>
24 #include <KokkosSparse_CrsMatrix.hpp>
29 template <
class Scalar,
32 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
34 :
public RowMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>,
35 public DistObject<char, LocalOrdinal, GlobalOrdinal, Node> {
55 const ArrayView<const GlobalOrdinal> &cols,
56 const ArrayView<const Scalar> &vals) = 0;
59 virtual void insertLocalValues(LocalOrdinal localRow,
const ArrayView<const LocalOrdinal> &cols,
const ArrayView<const Scalar> &vals) = 0;
62 virtual void replaceGlobalValues(GlobalOrdinal globalRow,
const ArrayView<const GlobalOrdinal> &cols,
const ArrayView<const Scalar> &vals) = 0;
65 virtual void replaceLocalValues(LocalOrdinal localRow,
const ArrayView<const LocalOrdinal> &cols,
const ArrayView<const Scalar> &vals) = 0;
71 virtual void scale(
const Scalar &alpha) = 0;
83 virtual void allocateAllValues(
size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) = 0;
86 virtual void setAllValues(
const ArrayRCP<size_t> &rowptr,
const ArrayRCP<LocalOrdinal> &colind,
const ArrayRCP<Scalar> &values) = 0;
89 virtual void getAllValues(ArrayRCP<const size_t> &rowptr, ArrayRCP<const LocalOrdinal> &colind, ArrayRCP<const Scalar> &values)
const = 0;
100 virtual void resumeFill(
const RCP<ParameterList> ¶ms = null) = 0;
106 virtual void fillComplete(
const RCP<ParameterList> ¶ms = null) = 0;
116 const RCP<ParameterList> ¶ms = Teuchos::null) = 0;
123 virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRowMap()
const = 0;
126 virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getColMap()
const = 0;
129 virtual RCP<const CrsGraph<LocalOrdinal, GlobalOrdinal, Node> >
getCrsGraph()
const = 0;
171 virtual typename ScalarTraits<Scalar>::magnitudeType
getFrobeniusNorm()
const = 0;
177 virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values)
const = 0;
180 virtual void getGlobalRowCopy(GlobalOrdinal GlobalRow,
const ArrayView<GlobalOrdinal> &indices,
const ArrayView<Scalar> &values,
size_t &numEntries)
const = 0;
183 virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values)
const = 0;
228 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;
247 virtual 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<
Import<LocalOrdinal, GlobalOrdinal, Node> > ®ionInterfaceImporter,
const Teuchos::ArrayRCP<LocalOrdinal> ®ionInterfaceLIDs)
const = 0;
250 virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getDomainMap()
const = 0;
253 virtual const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
getRangeMap()
const = 0;
264 virtual void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const = 0;
275 #ifdef HAVE_XPETRA_TPETRA
289 typename local_graph_type::size_type>;
294 virtual void setAllValues(
const typename local_matrix_type::row_map_type &ptr,
295 const typename local_graph_type::entries_type::non_const_type &ind,
296 const typename local_matrix_type::values_type &val) = 0;
299 #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."
311 virtual void getLocalRowCopy(LocalOrdinal LocalRow,
const ArrayView<LocalOrdinal> &Indices,
const ArrayView<Scalar> &Values,
size_t &NumEntries)
const = 0;
328 #define XPETRA_CRSMATRIX_SHORT
329 #endif // XPETRA_CRSMATRIX_HPP
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const =0
Returns the Frobenius norm of the matrix.
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 RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const =0
Returns the CrsGraph associated with this matrix.
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 void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)=0
Replace the diagonal entries of the matrix.
virtual const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const =0
Returns the Map that describes the row distribution in this 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 size_t getLocalNumCols() const =0
Returns the number of matrix columns owned on the calling node.
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 size_t getLocalNumEntries() const =0
Returns the local number of entries 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 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 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 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.
typename node_type::device_type execution_space
virtual global_size_t getGlobalNumCols() const =0
Number of global columns in the matrix.
Kokkos::StaticCrsGraph< int, Kokkos::LayoutLeft, execution_space, void, size_t > local_graph_type
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 size_t getLocalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
virtual local_matrix_type::HostMirror getLocalMatrixHost() const =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.
typename Kokkos::ArithTraits< double >::val_type impl_scalar_type
size_t global_size_t
Global size_t object.
virtual size_t getLocalNumRows() const =0
Returns the number of matrix rows owned on the calling node.
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 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 const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const =0
Returns the Map that describes the column distribution in this matrix.
virtual local_matrix_type getLocalMatrixDevice() const =0
virtual void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const =0
Get offsets of the diagonal entries in the matrix.
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 bool isGloballyIndexed() const =0
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
KokkosSparse::CrsMatrix< impl_scalar_type, int, 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...
virtual bool isFillActive() const =0
Returns true if the matrix is in edit mode.
virtual LocalOrdinal GetStorageBlockSize() const =0
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
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.