10#ifndef XPETRA_REORDEREDBLOCKEDMULTIVECTOR_HPP
11#define XPETRA_REORDEREDBLOCKEDMULTIVECTOR_HPP
13#include <Tpetra_KokkosCompat_DefaultNode.hpp>
21#include "Xpetra_BlockedMap.hpp"
22#include "Xpetra_BlockedMultiVector.hpp"
32template <
class Scalar,
35 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
44#undef XPETRA_REORDEREDBLOCKEDMULTIVECTOR_SHORT
82 size_t numBlocks = brm->GetNumBlocks();
91 map = bmap->getMap(Teuchos::as<size_t>(leaf->GetIndex()),
false);
94 std::vector<Teuchos::RCP<const Map>> subMaps(numBlocks, Teuchos::null);
96 for (
size_t i = 0; i < numBlocks; i++) {
113 std::string
description()
const {
return "ReorderedBlockedMultiVector"; }
129template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
137 size_t numBlocks = brm->GetNumBlocks();
141 if (numBlocks == 0) {
145 map = bmap->getMap(Teuchos::as<size_t>(leaf->GetIndex()), bThyraMode);
148 std::vector<Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>> subMaps(numBlocks, Teuchos::null);
150 for (
size_t i = 0; i < numBlocks; i++) {
174template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
184 size_t rowSz = rowMgr->GetNumBlocks();
199 if (subBVec == Teuchos::null) {
214 std::vector<Teuchos::RCP<const Map>> rowSubMaps(1, submap);
218 rbvec->setMultiVector(0, Teuchos::rcp_const_cast<MultiVector>(vec),
false);
229 std::vector<Teuchos::RCP<const Map>> rowSubMaps(rowSz, Teuchos::null);
230 for (
size_t i = 0; i < rowSz; i++) {
239 for (
size_t i = 0; i < rowSz; i++) {
242 rbvec->setMultiVector(i, Teuchos::rcp_const_cast<MultiVector>(subvec),
false);
248template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
259 size_t rowSz = rowMgr->GetNumBlocks();
274 if (bbvec == Teuchos::null) {
281 std::vector<Teuchos::RCP<const Map>> rowXpSubMaps(1, xpsubmap);
282 std::vector<Teuchos::RCP<const Map>> rowTySubMaps(1, thysubmap);
289 rbvec->setMultiVector(0, vec,
true);
292 rbvec = Teuchos::rcp_dynamic_cast<BlockedMultiVector>(vec);
298 std::vector<Teuchos::RCP<const Map>> rowXpSubMaps(rowSz, Teuchos::null);
299 std::vector<Teuchos::RCP<const Map>> rowTySubMaps(rowSz, Teuchos::null);
300 for (
size_t i = 0; i < rowSz; i++) {
313 for (
size_t i = 0; i < rowSz; i++) {
316 rbvec->setMultiVector(i, Teuchos::rcp_const_cast<MultiVector>(subvec),
true);
322template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 if (bvec->getBlockedMap()->getThyraMode() ==
false) {
334template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
339 if (bvec->getBlockedMap()->getThyraMode() ==
false) {
340 rbvec =
mergeSubBlocks(brm, Teuchos::rcp_const_cast<const BlockedMultiVector>(bvec));
345 return Teuchos::rcp_const_cast<MultiVector>(rbvec);
350#define XPETRA_REORDEREDBLOCKEDMULTIVECTOR_SHORT
static const EVerbosityLevel verbLevel_default
virtual size_t getNumVectors() const
Number of columns in the multivector.
Xpetra utility class for common map-related routines.
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.
std::string description() const
Return a simple one-line description of this object.
GlobalOrdinal global_ordinal_type
LocalOrdinal local_ordinal_type
virtual ~ReorderedBlockedMultiVector()
Destructor.
Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > fullVec_
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.
ReorderedBlockedMultiVector(Teuchos::RCP< const BlockedMap > &rangeMap, Teuchos::RCP< const Xpetra::BlockReorderManager > brm, Teuchos::RCP< const Xpetra::BlockedMultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > bvec)
Constructor.
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > mergeSubBlockMaps(Teuchos::RCP< const Xpetra::BlockReorderManager > brm)
Teuchos::RCP< const Xpetra::BlockReorderManager > brm_
#define TEUCHOS_ASSERT(assertion_test)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
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::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)
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)