12 #ifndef XPETRA_MATRIX_DECL_HPP
13 #define XPETRA_MATRIX_DECL_HPP
15 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
20 #include "Xpetra_MultiVector.hpp"
26 #include "Xpetra_StridedMap.hpp"
27 #include "Xpetra_StridedMapFactory.hpp"
29 #include <Teuchos_SerialDenseMatrix.hpp>
30 #include <Teuchos_Hashtable.hpp>
56 template <
class Scalar,
59 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
73 #ifdef HAVE_XPETRA_TPETRA
92 void CreateView(
const viewLabel_t viewLabel,
const RCP<const Matrix> &A,
bool transposeA =
false,
const RCP<const Matrix> &B = Teuchos::null,
bool transposeB =
false);
95 void PrintViews(Teuchos::FancyOStream &out)
const;
126 virtual void insertGlobalValues(GlobalOrdinal globalRow,
const ArrayView<const GlobalOrdinal> &cols,
const ArrayView<const Scalar> &vals) = 0;
136 virtual void insertLocalValues(LocalOrdinal localRow,
const ArrayView<const LocalOrdinal> &cols,
const ArrayView<const Scalar> &vals) = 0;
145 const ArrayView<const GlobalOrdinal> &cols,
146 const ArrayView<const Scalar> &vals) = 0;
153 const ArrayView<const LocalOrdinal> &cols,
154 const ArrayView<const Scalar> &vals) = 0;
160 virtual void scale(
const Scalar &alpha) = 0;
175 virtual void resumeFill(
const RCP<ParameterList> ¶ms = null) = 0;
188 virtual void fillComplete(
const RCP<const Map> &domainMap,
const RCP<const Map> &rangeMap,
const RCP<ParameterList> ¶ms = null) = 0;
204 virtual void fillComplete(
const RCP<ParameterList> ¶ms = null) = 0;
212 virtual const RCP<const Map> &
getRowMap()
const;
219 virtual const RCP<const Map> &
getColMap()
const;
284 const ArrayView<LocalOrdinal> &Indices,
285 const ArrayView<Scalar> &Values,
286 size_t &NumEntries)
const = 0;
298 virtual void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values)
const = 0;
310 virtual void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values)
const = 0;
318 virtual typename ScalarTraits<Scalar>::magnitudeType
getFrobeniusNorm()
const = 0;
335 virtual const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> >
getMap()
const = 0;
367 virtual void describe(Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default)
const = 0;
381 virtual RCP<const CrsGraph>
getCrsGraph()
const = 0;
389 Teuchos::ETransp mode,
392 bool sumInterfaceValues,
394 const Teuchos::ArrayRCP<LocalOrdinal> ®ionInterfaceLIDs)
const = 0;
425 #ifdef HAVE_XPETRA_TPETRA
430 #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."
450 #define XPETRA_MATRIX_SHORT
451 #endif // XPETRA_MATRIX_DECL_HPP
virtual void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using local IDs.
const viewLabel_t & GetCurrentViewLabel() const
virtual global_size_t getGlobalNumRows() const =0
Returns the number of global rows in this matrix.
virtual const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
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.
virtual ~Matrix()
Destructor.
void RemoveView(const viewLabel_t viewLabel)
virtual global_size_t getGlobalNumCols() const =0
Returns the number of global columns in the matrix.
virtual const RCP< const Map > & getRowMap() const
Returns the Map that describes the row distribution in this matrix.
LocalOrdinal local_ordinal_type
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const =0
Get Frobenius norm of the matrix.
virtual void resumeFill(const RCP< ParameterList > ¶ms=null)=0
virtual LocalOrdinal GetStorageBlockSize() const =0
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
virtual size_t getGlobalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on all nodes.
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
virtual const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const =0
Implements DistObject interface.
virtual bool isGloballyIndexed() const =0
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
const viewLabel_t & GetDefaultViewLabel() const
virtual void doExport(const Matrix &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0
Export.
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
Computes the matrix-multivector multiplication for region layout matrices.
virtual void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null)=0
Signal that data entry is complete, specifying domain and range maps.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
virtual void setObjectLabel(const std::string &objectLabel)=0
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 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 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...
viewLabel_t defaultViewLabel_
virtual void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Insert matrix entries, using global IDs.
void SetFixedBlockSize(LocalOrdinal blksize, GlobalOrdinal offset=0)
virtual bool haveGlobalConstants() const =0
Returns true if globalConstants have been computed; false otherwise.
virtual size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const =0
Returns the current number of entries in the specified global row.
virtual void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Right scale matrix using the given vector entries.
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.
bool IsView(const viewLabel_t viewLabel) const
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
virtual size_t getLocalNumEntries() const =0
Returns the local number of entries in this matrix.
virtual bool isFillComplete() const =0
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
virtual void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using local IDs.
viewLabel_t currentViewLabel_
virtual void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)=0
Left scale matrix using the given vector entries.
void PrintViews(Teuchos::FancyOStream &out) const
Print all of the views associated with the Matrix.
LocalOrdinal GetFixedBlockSize() const
virtual size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const =0
Returns the current number of entries on this node in the specified local row.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
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 void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)=0
Replace matrix entries, using global IDs.
virtual void doImport(const Matrix &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)=0
Import.
virtual RCP< const CrsGraph > getCrsGraph() const =0
Returns the CrsGraph associated with this matrix.
CrsMatrix::local_matrix_type local_matrix_type
virtual std::string description() const =0
Return a simple one-line description of this object.
virtual bool isLocallyIndexed() const =0
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
virtual Scalar GetMaxEigenvalueEstimate() const
virtual local_matrix_type::HostMirror getLocalMatrixHost() const =0
virtual void SetMaxEigenvalueEstimate(Scalar const &sigma)
bool IsFixedBlockSizeSet() const
Returns true, if SetFixedBlockSize has been called before.
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 global_size_t getGlobalNumEntries() const =0
Returns the global number of entries in this matrix.
virtual size_t getLocalMaxNumRowEntries() const =0
Returns the maximum number of entries across all rows/columns on this node.
const viewLabel_t SwitchToView(const viewLabel_t viewLabel)
CombineMode
Xpetra::Combine Mode enumerable type.
virtual void scale(const Scalar &alpha)=0
Scale the current values of a matrix, this = alpha*this.
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
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...
virtual void setAllToScalar(const Scalar &alpha)=0
Set all matrix entries equal to scalar.
const viewLabel_t SwitchToDefaultView()
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
virtual local_matrix_type getLocalMatrixDevice() const =0
Xpetra-specific matrix class.
virtual bool hasCrsGraph() const =0
Supports the getCrsGraph() call.
GlobalOrdinal global_ordinal_type