46 #ifndef XPETRA_REORDEREDBLOCKEDCRSMATRIX_HPP
47 #define XPETRA_REORDEREDBLOCKEDCRSMATRIX_HPP
69 template <
class Scalar,
82 #undef XPETRA_REORDEREDBLOCKEDCRSMATRIX_SHORT
132 map = fullRangeMapExtractor->getMap(Teuchos::as<size_t>(leaf->
GetIndex()),
false);
135 std::vector<Teuchos::RCP<const Map> > subMaps (numBlocks, Teuchos::null);
137 for(
size_t i = 0; i < numBlocks; i++) {
178 bool bCopyResultX =
false;
179 bool bCopyResultY =
false;
240 if (bCopyResultX ==
true) {
245 if (bCopyResultY ==
true) {
268 std::string
description()
const {
return "ReorderedBlockedCrsMatrix"; }
277 out <<
"ReorderedBlockMatrix is fillComplete" << std::endl;
279 out <<
"fullRowMap" << std::endl;
286 out <<
"Xpetra::ReorderedBlockedCrsMatrix is NOT fillComplete" << std::endl;
291 out <<
"Block(" << r <<
"," << c <<
")" << std::endl;
305 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
321 map = fullRangeMapExtractor->getMap(Teuchos::as<size_t>(leaf->
GetIndex()), bThyraMode);
324 std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > > subMaps (numBlocks, Teuchos::null);
326 for(
size_t i = 0; i < numBlocks; i++) {
351 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
367 if(rowSz == 0 && colSz == 0) {
375 if (mat == Teuchos::null)
return Teuchos::null;
379 if(matwrap != Teuchos::null) {
384 std::vector<Teuchos::RCP<const Map> > rowSubMaps (1, submap);
389 std::vector<Teuchos::RCP<const Map> > colSubMaps (1, submap2);
392 rbmat =
Teuchos::rcp(
new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
396 rbmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
405 std::vector<Teuchos::RCP<const Map> > rowSubMaps (rowSz, Teuchos::null);
406 for(
size_t i = 0; i < rowSz; i++) {
412 rgMapExtractor =
Teuchos::rcp(
new MapExtractor(rgMergedSubMaps, rowSubMaps,
false));
419 std::vector<Teuchos::RCP<const Map> > rowSubMaps (1, submap);
420 rgMapExtractor =
Teuchos::rcp(
new MapExtractor(submap, rowSubMaps,
false));
425 std::vector<Teuchos::RCP<const Map> > colSubMaps (colSz, Teuchos::null);
426 for(
size_t j = 0; j < colSz; j++) {
432 doMapExtractor =
Teuchos::rcp(
new MapExtractor(doMergedSubMaps, colSubMaps,
false));
439 std::vector<Teuchos::RCP<const Map> > colSubMaps (1, submap);
440 doMapExtractor =
Teuchos::rcp(
new MapExtractor(submap, colSubMaps,
false));
443 rbmat =
Teuchos::rcp(
new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
447 if (rowSz == 0 && colSz > 0) {
448 for(
size_t j = 0; j < colSz; j++) {
451 rbmat->
setMatrix(0,j,Teuchos::rcp_const_cast<Matrix>(submat));
454 }
else if (rowSz > 0 && colSz == 0) {
455 for(
size_t i = 0; i < rowSz; i++) {
458 rbmat->
setMatrix(i,0,Teuchos::rcp_const_cast<Matrix>(submat));
462 for(
size_t i = 0; i < rowSz; i++) {
464 for(
size_t j = 0; j < colSz; j++) {
467 rbmat->
setMatrix(i,j,Teuchos::rcp_const_cast<Matrix>(submat));
480 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
489 TEUCHOS_ASSERT(bmat->getRangeMapExtractor()->getThyraMode() ==
true);
490 TEUCHOS_ASSERT(bmat->getDomainMapExtractor()->getThyraMode() ==
true);
498 if(rowSz == 0 && colSz == 0) {
506 if(mat == Teuchos::null)
return Teuchos::null;
510 if(matwrap != Teuchos::null) {
517 std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (1, xpsubmap);
518 std::vector<Teuchos::RCP<const Map> > rowTySubMaps (1, thysubmap);
526 std::vector<Teuchos::RCP<const Map> > colXpSubMaps (1, xpsubmap2);
527 std::vector<Teuchos::RCP<const Map> > colTySubMaps (1, tysubmap2);
533 rbmat =
Teuchos::rcp(
new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
537 rbmat = Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(mat);
546 std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (rowSz, Teuchos::null);
547 std::vector<Teuchos::RCP<const Map> > rowTySubMaps (rowSz, Teuchos::null);
548 for(
size_t i = 0; i < rowSz; i++) {
557 rgMapExtractor =
Teuchos::rcp(
new MapExtractor(rowXpSubMaps, rowTySubMaps));
564 std::vector<Teuchos::RCP<const Map> > rowXpSubMaps (1, xpsubmap);
565 std::vector<Teuchos::RCP<const Map> > rowTySubMaps (1, thysubmap);
567 rgMapExtractor =
Teuchos::rcp(
new MapExtractor(rowXpSubMaps, rowTySubMaps));
572 std::vector<Teuchos::RCP<const Map> > colXpSubMaps (colSz, Teuchos::null);
573 std::vector<Teuchos::RCP<const Map> > colTySubMaps (colSz, Teuchos::null);
574 for(
size_t j = 0; j < colSz; j++) {
583 doMapExtractor =
Teuchos::rcp(
new MapExtractor(colXpSubMaps,colTySubMaps));
590 std::vector<Teuchos::RCP<const Map> > colXpSubMaps (1, xpsubmap);
591 std::vector<Teuchos::RCP<const Map> > colTySubMaps (1, tysubmap);
593 doMapExtractor =
Teuchos::rcp(
new MapExtractor(colXpSubMaps, colTySubMaps));
597 rbmat =
Teuchos::rcp(
new ReorderedBlockedCrsMatrix(rgMapExtractor,doMapExtractor, 33, rowMgr,bmat));
601 if (rowSz == 0 && colSz > 0) {
602 for(
size_t j = 0; j < colSz; j++) {
605 rbmat->
setMatrix(0,j,Teuchos::rcp_const_cast<Matrix>(submat));
608 }
else if (rowSz > 0 && colSz == 0) {
609 for(
size_t i = 0; i < rowSz; i++) {
612 rbmat->
setMatrix(i,0,Teuchos::rcp_const_cast<Matrix>(submat));
616 for(
size_t i = 0; i < rowSz; i++) {
618 for(
size_t j = 0; j < colSz; j++) {
621 rbmat->
setMatrix(i,j,Teuchos::rcp_const_cast<Matrix>(submat));
633 template<
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
635 TEUCHOS_ASSERT(bmat->getRangeMapExtractor()->getThyraMode() == bmat->getDomainMapExtractor()->getThyraMode());
637 if(bmat->getRangeMapExtractor()->getThyraMode() ==
false) {
649 #define XPETRA_REORDEREDBLOCKEDCRSMATRIX_SHORT
virtual size_t getNodeNumEntries() const =0
Returns the local number of entries in this matrix.
bool is_null(const boost::shared_ptr< T > &p)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat, bool bThyraMode)
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete.
RCP< const Map > getRangeMap() const
Returns the Map associated with the range of this operator.
Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > fullOp_
void swap(RCP< T > &r_ptr)
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
Xpetra utility class for common map-related routines.
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.
Teuchos::RCP< const Xpetra::BlockReorderManager > getBlockReorderManager()
Returns internal BlockReorderManager object.
static Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > concatenateMaps(const std::vector< Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > > &subMaps)
Helper function to concatenate several maps.
GlobalOrdinal global_ordinal_type
BlockedCrsMatrix(const Teuchos::RCP< const BlockedMap > &rangeMaps, const Teuchos::RCP< const BlockedMap > &domainMaps, size_t npr, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
int GetIndex() const
Get the index that is stored in this block/leaf.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlocksThyra(Teuchos::RCP< const Xpetra::BlockReorderManager > rowMgr, Teuchos::RCP< const Xpetra::BlockReorderManager > colMgr, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat)
virtual void update(const Scalar &alpha, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, const Scalar &beta)=0
Update multi-vector values with scaled values of A, this = beta*this + alpha*A.
virtual size_t GetNumBlocks() const
Returns number of subblocks.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlocks(Teuchos::RCP< const Xpetra::BlockReorderManager > rowMgr, Teuchos::RCP< const Xpetra::BlockReorderManager > colMgr, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< Matrix > getMatrix(size_t r, size_t c) const
return block (r,c)
virtual const Teuchos::RCP< BlockReorderManager > GetBlock(int blockIndex)
Get a particular block. If there is no block at this index location return a new one.
virtual size_t Cols() const
number of column blocks
static const EVerbosityLevel verbLevel_default
std::string description() const
Return a simple one-line description of this object.
Teuchos::RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > buildReorderedBlockedCrsMatrix(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bmat)
virtual ~ReorderedBlockedCrsMatrix()
Destructor.
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm)
void setMatrix(size_t r, size_t c, Teuchos::RCP< Matrix > mat)
set matrix block
LocalOrdinal local_ordinal_type
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.
Concrete implementation of Xpetra::Matrix.
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, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor.
virtual size_t Rows() const
number of row blocks
Teuchos::RCP< const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > buildReorderedBlockedMultiVector(Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bvec)
#define TEUCHOS_ASSERT(assertion_test)
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
Xpetra-specific matrix class.
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.
Teuchos::RCP< const Xpetra::BlockReorderManager > brm_