10 #include "Teko_TpetraBlockedMappingStrategy.hpp"
11 #include "Teko_TpetraHelpers.hpp"
13 #include "Thyra_TpetraThyraWrappers.hpp"
14 #include "Thyra_TpetraLinearOp.hpp"
15 #include "Thyra_DefaultProductMultiVector.hpp"
16 #include "Thyra_DefaultProductVectorSpace.hpp"
17 #include "Thyra_DefaultSpmdMultiVector.hpp"
18 #include "Thyra_DefaultBlockedLinearOp.hpp"
22 using Teuchos::rcp_dynamic_cast;
25 namespace TpetraHelpers {
37 TpetraBlockedMappingStrategy::TpetraBlockedMappingStrategy(
38 const std::vector<std::vector<GO> >& vars,
39 const Teuchos::RCP<
const Tpetra::Map<LO, GO, NT> >& map,
const Teuchos::Comm<int>& comm) {
42 buildBlockTransferData(vars, rangeMap_, comm);
53 void TpetraBlockedMappingStrategy::copyTpetraIntoThyra(
54 const Tpetra::MultiVector<ST, LO, GO, NT>& X,
55 const Teuchos::Ptr<Thyra::MultiVectorBase<ST> >& thyra_X)
const {
56 int count = X.getNumVectors();
58 std::vector<RCP<Tpetra::MultiVector<ST, LO, GO, NT> > > subX;
61 Blocking::buildSubVectors(blockMaps_, subX, count);
64 Blocking::one2many(subX, X, blockImport_);
67 Teuchos::Array<RCP<Thyra::MultiVectorBase<ST> > > thyra_subX;
68 Teuchos::Ptr<Thyra::ProductMultiVectorBase<ST> > prod_X =
69 Teuchos::ptr_dynamic_cast<Thyra::ProductMultiVectorBase<ST> >(thyra_X);
70 for (
unsigned int i = 0; i < blockMaps_.size(); i++) {
71 RCP<Thyra::TpetraMultiVector<ST, LO, GO, NT> > vec =
72 rcp_dynamic_cast<Thyra::TpetraMultiVector<ST, LO, GO, NT> >(
73 prod_X->getNonconstMultiVectorBlock(i),
true);
75 fillDefaultSpmdMultiVector(vec, subX[i]);
87 void TpetraBlockedMappingStrategy::copyThyraIntoTpetra(
88 const RCP<
const Thyra::MultiVectorBase<ST> >& thyra_Y,
89 Tpetra::MultiVector<ST, LO, GO, NT>& Y)
const {
90 std::vector<RCP<const Tpetra::MultiVector<ST, LO, GO, NT> > > subY;
91 RCP<const Thyra::DefaultProductMultiVector<ST> > prod_Y =
92 rcp_dynamic_cast<
const Thyra::DefaultProductMultiVector<ST> >(thyra_Y);
95 for (
unsigned int i = 0; i < blockMaps_.size(); i++) {
96 RCP<const Thyra::TpetraMultiVector<ST, LO, GO, NT> > tmv =
97 rcp_dynamic_cast<
const Thyra::TpetraMultiVector<ST, LO, GO, NT> >(
98 prod_Y->getMultiVectorBlock(i),
true);
99 subY.push_back(tmv->getConstTpetraMultiVector());
106 Blocking::many2one(Y, subY, blockExport_);
121 void TpetraBlockedMappingStrategy::buildBlockTransferData(
122 const std::vector<std::vector<GO> >& vars,
123 const Teuchos::RCP<
const Tpetra::Map<LO, GO, NT> >& baseMap,
const Teuchos::Comm<int>& comm) {
125 for (std::size_t i = 0; i < vars.size(); i++) {
127 Blocking::MapPair mapPair = Blocking::buildSubMap(vars[i], comm);
128 Blocking::ImExPair iePair = Blocking::buildExportImport(*baseMap, mapPair);
130 blockMaps_.push_back(mapPair);
131 blockImport_.push_back(iePair.first);
132 blockExport_.push_back(iePair.second);
146 const Teuchos::RCP<Thyra::BlockedLinearOpBase<ST> >
147 TpetraBlockedMappingStrategy::buildBlockedThyraOp(
148 const RCP<
const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
149 const std::string& label)
const {
150 int dim = blockMaps_.size();
152 RCP<Thyra::DefaultBlockedLinearOp<ST> > A = Thyra::defaultBlockedLinearOp<ST>();
154 A->beginBlockFill(dim, dim);
155 for (
int i = 0; i < dim; i++) {
156 for (
int j = 0; j < dim; j++) {
158 std::stringstream ss;
159 ss << label <<
"_" << i <<
"," << j;
162 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > blk =
163 Blocking::buildSubBlock(i, j, crsContent, blockMaps_);
164 A->setNonconstBlock(i, j,
165 Thyra::tpetraLinearOp<ST, LO, GO, NT>(
166 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(blk->getRangeMap()),
167 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(blk->getDomainMap()), blk));
185 void TpetraBlockedMappingStrategy::rebuildBlockedThyraOp(
186 const RCP<
const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
187 const RCP<Thyra::BlockedLinearOpBase<ST> >& A)
const {
188 int dim = blockMaps_.size();
190 for (
int i = 0; i < dim; i++) {
191 for (
int j = 0; j < dim; j++) {
193 RCP<Thyra::LinearOpBase<ST> > Aij = A->getNonconstBlock(i, j);
194 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT> > tAij =
195 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT> >(Aij,
true);
196 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > eAij =
197 rcp_dynamic_cast<Tpetra::CrsMatrix<ST, LO, GO, NT> >(tAij->getTpetraOperator(),
true);
200 Blocking::rebuildSubBlock(i, j, crsContent, blockMaps_, *eAij);