Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

#include <Xpetra_ReorderedBlockedCrsMatrix.hpp>

Inheritance diagram for Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Public Types

typedef Scalar scalar_type
 
typedef LocalOrdinal local_ordinal_type
 
typedef GlobalOrdinal global_ordinal_type
 
typedef Node node_type
 
- Public Types inherited from Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
typedef Scalar scalar_type
 
typedef LocalOrdinal local_ordinal_type
 
typedef GlobalOrdinal global_ordinal_type
 
typedef Node node_type
 
typedef
CrsMatrix::local_matrix_type 
local_matrix_type
 
- Public Types inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
typedef Scalar scalar_type
 
typedef LocalOrdinal local_ordinal_type
 
typedef GlobalOrdinal global_ordinal_type
 
typedef Node node_type
 
typedef
CrsMatrix::local_matrix_type 
local_matrix_type
 
- Public Types inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >
typedef Scalar scalar_type
 The type of the entries of the input and output multivectors. More...
 
typedef LocalOrdinal local_ordinal_type
 The local index type. More...
 
typedef GlobalOrdinal global_ordinal_type
 The global index type. More...
 
typedef Node node_type
 The Kokkos Node type. More...
 

Private Member Functions

Teuchos::RCP< const
Xpetra::Map< LocalOrdinal,
GlobalOrdinal, Node > > 
mergeSubBlockMaps (Teuchos::RCP< const Xpetra::BlockReorderManager > brm)
 

Private Attributes

Teuchos::RCP< const
Xpetra::BlockReorderManager
brm_
 
Teuchos::RCP< const
Xpetra::BlockedCrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
fullOp_
 

Constructor/Destructor Methods

 ReorderedBlockedCrsMatrix (Teuchos::RCP< const MapExtractor > &rangeMaps, Teuchos::RCP< const MapExtractor > &domainMaps, size_t npr, Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> bmat)
 Constructor. More...
 
virtual ~ReorderedBlockedCrsMatrix ()
 Destructor. More...
 

Methods implementing Matrix

virtual void apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node >> &regionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > &regionInterfaceLIDs) const
 sparse matrix-multivector multiplication for the region layout matrices (currently no blocked implementation) More...
 
virtual 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. More...
 

Access functions

Teuchos::RCP< const
Xpetra::BlockReorderManager
getBlockReorderManager ()
 Returns internal BlockReorderManager object. More...
 
Teuchos::RCP< const
Xpetra::BlockedCrsMatrix
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
getBlockedCrsMatrix ()
 Returns internal unmodified BlockedCrsMatrix object. More...
 

Overridden from Teuchos::Describable

std::string description () const
 Return a simple one-line description of this object. More...
 
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. More...
 

Additional Inherited Members

- Public Member Functions inherited from Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
global_size_t getGlobalNumRows () const
 Returns the number of global rows. More...
 
global_size_t getGlobalNumCols () const
 Returns the number of global columns in the matrix. More...
 
size_t getLocalNumRows () const
 Returns the number of matrix rows owned on the calling node. More...
 
global_size_t getGlobalNumEntries () const
 Returns the global number of entries in this matrix. More...
 
size_t getLocalNumEntries () const
 Returns the local number of entries in this matrix. More...
 
size_t getNumEntriesInLocalRow (LocalOrdinal localRow) const
 Returns the current number of entries on this node in the specified local row. More...
 
size_t getNumEntriesInGlobalRow (GlobalOrdinal globalRow) const
 Returns the current number of entries in the specified (locally owned) global row. More...
 
size_t getGlobalMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on all nodes. More...
 
size_t getLocalMaxNumRowEntries () const
 Returns the maximum number of entries across all rows/columns on this node. More...
 
bool isLocallyIndexed () const
 If matrix indices of all matrix blocks are in the local range, this function returns true. Otherwise, this function returns false. More...
 
bool isGloballyIndexed () const
 If matrix indices are in the global range, this function returns true. Otherwise, this function returns false. More...
 
bool isFillComplete () const
 Returns true if fillComplete() has been called and the matrix is in compute mode. More...
 
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 calling routine. More...
 
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. More...
 
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. More...
 
void getLocalDiagCopy (Vector &diag) const
 Get a copy of the diagonal entries owned by this node, with local row indices. More...
 
void leftScale (const Vector &x)
 Left scale matrix using the given vector entries. More...
 
void rightScale (const Vector &x)
 Right scale matrix using the given vector entries. More...
 
virtual ScalarTraits< Scalar >
::magnitudeType 
getFrobeniusNorm () const
 Get Frobenius norm of the matrix. More...
 
virtual bool haveGlobalConstants () const
 Returns true if globalConstants have been computed; false otherwise. More...
 
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) More...
 
const Teuchos::RCP< const MapgetMap () const
 Implements DistObject interface. More...
 
void doImport (const Matrix &source, const Import &importer, CombineMode CM)
 Import. More...
 
void doExport (const Matrix &dest, const Import &importer, CombineMode CM)
 Export. More...
 
void doImport (const Matrix &source, const Export &exporter, CombineMode CM)
 Import (using an Exporter). More...
 
void doExport (const Matrix &dest, const Export &exporter, CombineMode CM)
 Export (using an Importer). More...
 
bool hasCrsGraph () const
 Supports the getCrsGraph() call. More...
 
RCP< const CrsGraphgetCrsGraph () const
 Returns the CrsGraph associated with this matrix. More...
 
local_matrix_type getLocalMatrixDevice () const
 Access the underlying local Kokkos::CrsMatrix object. More...
 
local_matrix_type::HostMirror getLocalMatrixHost () const
 Access the underlying local Kokkos::CrsMatrix object. More...
 
LocalOrdinal GetStorageBlockSize () const
 Returns the block size of the storage mechanism. More...
 
void residual (const MultiVector &X, const MultiVector &B, MultiVector &R) const
 Compute a residual R = B - (*this) * X. More...
 
 BlockedCrsMatrix (const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t numEntriesPerRow)
 Constructor. More...
 
 BlockedCrsMatrix (Teuchos::RCP< const MapExtractor > &rangeMapExtractor, Teuchos::RCP< const MapExtractor > &domainMapExtractor, size_t numEntriesPerRow)
 Constructor. More...
 
virtual ~BlockedCrsMatrix ()
 Destructor. More...
 
void insertGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
 Insert matrix entries, using global IDs. More...
 
void insertLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
 Insert matrix entries, using local IDs. More...
 
void removeEmptyProcessesInPlace (const Teuchos::RCP< const Map > &newMap)
 
void replaceGlobalValues (GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
 Replace matrix entries, using global IDs. More...
 
void replaceLocalValues (LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
 Replace matrix entries, using local IDs. More...
 
virtual void setAllToScalar (const Scalar &alpha)
 Set all matrix entries equal to scalar. More...
 
void scale (const Scalar &alpha)
 Scale the current values of a matrix, this = alpha*this. More...
 
void resumeFill (const RCP< ParameterList > &params=null)
 
void fillComplete (const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > &params=null)
 Signal that data entry is complete. More...
 
void fillComplete (const RCP< ParameterList > &params=null)
 Signal that data entry is complete. More...
 
virtual void apply (const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > > &regionInterfaceImporter, const Teuchos::ArrayRCP< LocalOrdinal > &regionInterfaceLIDs) const
 sparse matrix-multivector multiplication for the region layout matrices (currently no blocked implementation) More...
 
RCP< const MapgetFullDomainMap () const
 Returns the Map associated with the full domain of this operator. More...
 
RCP< const BlockedMapgetBlockedDomainMap () const
 Returns the BlockedMap associated with the domain of this operator. More...
 
const RCP< const MapgetDomainMap () const
 Returns the Map associated with the domain of this operator. More...
 
RCP< const MapgetDomainMap (size_t i) const
 Returns the Map associated with the i'th block domain of this operator. More...
 
RCP< const MapgetDomainMap (size_t i, bool bThyraMode) const
 Returns the Map associated with the i'th block domain of this operator. More...
 
RCP< const MapgetFullRangeMap () const
 Returns the Map associated with the full range of this operator. More...
 
RCP< const BlockedMapgetBlockedRangeMap () const
 Returns the BlockedMap associated with the range of this operator. More...
 
const RCP< const MapgetRangeMap () const
 Returns the Map associated with the range of this operator. More...
 
RCP< const MapgetRangeMap (size_t i) const
 Returns the Map associated with the i'th block range of this operator. More...
 
RCP< const MapgetRangeMap (size_t i, bool bThyraMode) const
 Returns the Map associated with the i'th block range of this operator. More...
 
RCP< const MapExtractorgetRangeMapExtractor () const
 Returns map extractor class for range map. More...
 
RCP< const MapExtractorgetDomainMapExtractor () const
 Returns map extractor for domain map. More...
 
void setObjectLabel (const std::string &objectLabel)
 
virtual bool isDiagonal () const
 
virtual size_t Rows () const
 number of row blocks More...
 
virtual size_t Cols () const
 number of column blocks More...
 
Teuchos::RCP< MatrixgetCrsMatrix () const
 return unwrap 1x1 blocked operators More...
 
Teuchos::RCP< MatrixgetInnermostCrsMatrix ()
 helper routine recursively returns the first inner-most non-null matrix block from a (nested) blocked operator More...
 
Teuchos::RCP< MatrixgetMatrix (size_t r, size_t c) const
 return block (r,c) More...
 
void setMatrix (size_t r, size_t c, Teuchos::RCP< Matrix > mat)
 set matrix block More...
 
Teuchos::RCP< MatrixMerge () const
 merge BlockedCrsMatrix blocks in a CrsMatrix More...
 
- Public Member Functions inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
void SetFixedBlockSize (LocalOrdinal blksize, GlobalOrdinal offset=0)
 
LocalOrdinal GetFixedBlockSize () const
 
bool IsFixedBlockSizeSet () const
 Returns true, if SetFixedBlockSize has been called before. More...
 
virtual void SetMaxEigenvalueEstimate (Scalar const &sigma)
 
virtual Scalar GetMaxEigenvalueEstimate () const
 
 Matrix ()
 
virtual ~Matrix ()
 Destructor. More...
 
void CreateView (viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
 
void CreateView (const viewLabel_t viewLabel, const RCP< const Matrix > &A, bool transposeA=false, const RCP< const Matrix > &B=Teuchos::null, bool transposeB=false)
 
void PrintViews (Teuchos::FancyOStream &out) const
 Print all of the views associated with the Matrix. More...
 
void RemoveView (const viewLabel_t viewLabel)
 
const viewLabel_t SwitchToView (const viewLabel_t viewLabel)
 
bool IsView (const viewLabel_t viewLabel) const
 
const viewLabel_t SwitchToDefaultView ()
 
const viewLabel_tGetDefaultViewLabel () const
 
const viewLabel_tGetCurrentViewLabel () const
 
virtual const RCP< const Map > & getRowMap () const
 Returns the Map that describes the row distribution in this matrix. More...
 
virtual const RCP< const Map > & getRowMap (viewLabel_t viewLabel) const
 Returns the Map that describes the row distribution in this matrix. More...
 
virtual const RCP< const Map > & getColMap () const
 Returns the Map that describes the column distribution in this matrix. This might be null until fillComplete() is called. More...
 
virtual const RCP< const Map > & getColMap (viewLabel_t viewLabel) const
 Returns the Map that describes the column distribution in this matrix. More...
 
- Public Member Functions inherited from Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual ~Operator ()
 
virtual void removeEmptyProcessesInPlace (const RCP< const map_type > &)
 
virtual bool hasTransposeApply () const
 Whether this operator supports applying the transpose or conjugate transpose. More...
 
- Protected Attributes inherited from Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >
Teuchos::Hashtable
< viewLabel_t, RCP< MatrixView > > 
operatorViewTable_
 
viewLabel_t defaultViewLabel_
 
viewLabel_t currentViewLabel_
 

Detailed Description

template<class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
class Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Definition at line 73 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

Member Typedef Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
typedef Scalar Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalar_type

Definition at line 75 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
typedef LocalOrdinal Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_ordinal_type

Definition at line 76 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
typedef GlobalOrdinal Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::global_ordinal_type

Definition at line 77 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
typedef Node Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::node_type

Definition at line 78 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

Constructor & Destructor Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ReorderedBlockedCrsMatrix ( Teuchos::RCP< const MapExtractor > &  rangeMaps,
Teuchos::RCP< const MapExtractor > &  domainMaps,
size_t  npr,
Teuchos::RCP< const Xpetra::BlockReorderManager brm,
Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >>  bmat 
)
inline

Constructor.

Parameters
rangeMapsrange maps for all blocks
domainMapsdomain maps for all blocks
nprextimated number of entries per row in each block(!)
brmof type BlockReorderManager
bmatoriginal full blocked operator (we keep the RCP to make sure all subblocks are available)

Definition at line 96 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
virtual Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~ReorderedBlockedCrsMatrix ( )
inlinevirtual

Destructor.

Definition at line 109 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::mergeSubBlockMaps ( Teuchos::RCP< const Xpetra::BlockReorderManager brm)
inlineprivate

Definition at line 114 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
virtual void Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const MultiVector X,
MultiVector Y,
Teuchos::ETransp  mode,
Scalar  alpha,
Scalar  beta,
bool  sumInterfaceValues,
const RCP< Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node >> &  regionInterfaceImporter,
const Teuchos::ArrayRCP< LocalOrdinal > &  regionInterfaceLIDs 
) const
inlinevirtual

sparse matrix-multivector multiplication for the region layout matrices (currently no blocked implementation)

Definition at line 149 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
virtual void Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::apply ( const MultiVector X,
MultiVector Y,
Teuchos::ETransp  mode = Teuchos::NO_TRANS,
Scalar  alpha = ScalarTraits<Scalar>::one(),
Scalar  beta = ScalarTraits<Scalar>::zero() 
) const
inlinevirtual

Computes the sparse matrix-multivector multiplication.

Performs \(Y = \alpha A^{\textrm{mode}} X + \beta Y\), with one special exceptions:

  • if beta == 0, apply() overwrites Y, so that any values in Y (including NaNs) are ignored.

Reimplemented from Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 156 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Teuchos::RCP<const Xpetra::BlockReorderManager> Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockReorderManager ( )
inline

Returns internal BlockReorderManager object.

Definition at line 253 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Teuchos::RCP<const Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::getBlockedCrsMatrix ( )
inline

Returns internal unmodified BlockedCrsMatrix object.

Definition at line 256 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
std::string Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::description ( ) const
inlinevirtual

Return a simple one-line description of this object.

Reimplemented from Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 264 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
void Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel = Teuchos::Describable::verbLevel_default 
) const
inlinevirtual

Print the object with some verbosity level to an FancyOStream object.

Reimplemented from Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 267 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

Member Data Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Teuchos::RCP<const Xpetra::BlockReorderManager> Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::brm_
private

Definition at line 293 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
Teuchos::RCP<const Xpetra::BlockedCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> > Xpetra::ReorderedBlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::fullOp_
private

Definition at line 294 of file Xpetra_ReorderedBlockedCrsMatrix.hpp.


The documentation for this class was generated from the following file: