10#ifndef XPETRA_BLOCKEDCRSMATRIX_DECL_HPP
11#define XPETRA_BLOCKEDCRSMATRIX_DECL_HPP
13#include <Tpetra_KokkosCompat_DefaultNode.hpp>
21#include "Xpetra_MapFactory.hpp"
22#include "Xpetra_MultiVector.hpp"
23#include "Xpetra_BlockedMultiVector.hpp"
24#include "Xpetra_MultiVectorFactory.hpp"
25#include "Xpetra_BlockedVector.hpp"
30#include "Xpetra_MapExtractor.hpp"
33#include "Xpetra_Matrix.hpp"
34#include "Xpetra_MatrixFactory.hpp"
35#include "Xpetra_CrsMatrixWrap.hpp"
37#ifdef HAVE_XPETRA_THYRA
38#include <Thyra_ProductVectorSpaceBase.hpp>
39#include <Thyra_VectorSpaceBase.hpp>
40#include <Thyra_LinearOpBase.hpp>
41#include <Thyra_BlockedLinearOpBase.hpp>
42#include <Thyra_PhysicallyBlockedLinearOpBase.hpp>
43#include "Xpetra_ThyraUtils.hpp"
46#include "Xpetra_VectorFactory.hpp"
54#ifdef HAVE_XPETRA_THYRA
55template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
61template <
class Scalar,
64 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
73#undef XPETRA_BLOCKEDCRSMATRIX_SHORT
88 size_t numEntriesPerRow);
100 size_t numEntriesPerRow);
102#ifdef HAVE_XPETRA_THYRA
202 void scale(
const Scalar& alpha);
313 size_t& NumEntries)
const;
509 virtual size_t Rows()
const;
512 virtual size_t Cols()
const;
537#if KOKKOS_VERSION >= 40799
543#ifdef HAVE_XPETRA_THYRA
569 void Add(
const Matrix& A,
const Scalar scalarA,
Matrix& B,
const Scalar scalarB)
const;
583#ifdef HAVE_XPETRA_THYRA
592#define XPETRA_BLOCKEDCRSMATRIX_SHORT
static const EVerbosityLevel verbLevel_default
Teuchos::RCP< const MapExtractor > domainmaps_
full domain map together with all partial domain maps
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
virtual ~BlockedCrsMatrix()
Destructor.
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 size_t Cols() const
number of column blocks
bool isLocallyIndexed() const
If matrix indices of all matrix blocks are in the local range, this function returns true....
Teuchos::RCP< const MapExtractor > rangemaps_
full range map together with all partial domain maps
const RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete.
bool bDomainThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
const Teuchos::RCP< const Map > getMap() const
Implements DistObject interface.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
std::vector< Teuchos::RCP< Matrix > > blocks_
row major matrix block storage
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
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.
RCP< const BlockedMap > getBlockedDomainMap() const
Returns the BlockedMap associated with the domain of this operator.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
GlobalOrdinal global_ordinal_type
Teuchos::RCP< Matrix > Merge() const
merge BlockedCrsMatrix blocks in a CrsMatrix
local_matrix_type::HostMirror getLocalMatrixHost() const
Access the underlying local Kokkos::CrsMatrix object.
void leftScale(const Vector &x)
Left scale matrix using the given vector entries.
const RCP< const Map > getDomainMap() const
Returns the Map associated with the domain of this operator.
Teuchos::RCP< Matrix > getCrsMatrix() const
return unwrap 1x1 blocked operators
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
void doImport(const Matrix &source, const Import &importer, CombineMode CM)
Import.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
void residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
bool is_diagonal_
If we're diagonal, a bunch of the extraction stuff should work.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
virtual bool isDiagonal() const
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
LocalOrdinal local_ordinal_type
void Add(const Matrix &A, const Scalar scalarA, Matrix &B, const Scalar scalarB) const
Add a Xpetra::CrsMatrix to another: B = B*scalarB + A*scalarA.
bool bRangeThyraMode_
boolean flag, which is true, if BlockedCrsMatrix has been created using Thyra-style numbering for sub...
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
virtual size_t Rows() const
number of row blocks
RCP< const MapExtractor > getRangeMapExtractor() const
Returns map extractor class for range map.
RCP< const BlockedMap > getBlockedRangeMap() const
Returns the BlockedMap associated with the range of this operator.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism.
virtual void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > ®ionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > ®ionInterfaceLIDs) const
sparse matrix-multivector multiplication for the region layout matrices (currently no blocked impleme...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified (locally owned) global row.
CrsMatrix::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.
virtual ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
virtual void bgs_apply(const MultiVector &X, MultiVector &Y, size_t row, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Special multiplication routine (for BGS/Jacobi smoother)
RCP< const Map > getFullDomainMap() const
Returns the Map associated with the full domain of this operator.
void rightScale(const Vector &x)
Right scale matrix using the given vector entries.
void setObjectLabel(const std::string &objectLabel)
RCP< const Map > getFullRangeMap() const
Returns the Map associated with the full range of this operator.
void resumeFill(const RCP< ParameterList > ¶ms=null)
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Teuchos::RCP< Matrix > getInnermostCrsMatrix()
helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked...
global_size_t getGlobalNumRows() const
Returns the number of global rows.
RCP< const MapExtractor > getDomainMapExtractor() const
Returns map extractor for domain map.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
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 removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
virtual 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...
void doExport(const Matrix &dest, const Import &importer, CombineMode CM)
Export.
std::string description() const
Return a simple one-line description of this object.
local_matrix_type getLocalMatrixDevice() const
Access the underlying local Kokkos::CrsMatrix object.
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...
Xpetra-specific matrix class.
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.