Teko Version of the Day
Loading...
Searching...
No Matches
Teko_TpetraBlockedMappingStrategy.cpp
1// @HEADER
2// *****************************************************************************
3// Teko: A package for block and physics based preconditioning
4//
5// Copyright 2010 NTESS and the Teko contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#include "Teko_TpetraBlockedMappingStrategy.hpp"
11#include "Teko_TpetraHelpers.hpp"
12
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"
19
20using Teuchos::RCP;
21using Teuchos::rcp;
22using Teuchos::rcp_dynamic_cast;
23
24namespace Teko {
25namespace TpetraHelpers {
26
27// Creates a strided mapping strategy. This class is useful
28// for breaking up nodally ordered matrices (i.e. the unknowns
29// in a FEM problem are ordered [u0,v0,p0,u1,v1,p1,...]). Current
30// implimentation only supports a fixed number of variables
31//
32// arguments:
33// vars - Number of different variables
34// map - original Tpetra::Map<LO,GO,NT> to be broken up
35// comm - Teuchos::RCP<Teuchos::Comm<int> > object related to the map
36//
37TpetraBlockedMappingStrategy::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) {
40 rangeMap_ = map;
41 domainMap_ = map;
42 buildBlockTransferData(vars, rangeMap_, comm);
43}
44
45// Virtual function defined in MappingStrategy. This copies
46// an Tpetra::MultiVector<ST,LO,GO,NT> into a Thyra::MultiVectorBase with
47// blocking handled by the strides defined in the constructor.
48//
49// arguments:
50// X - source Tpetra::MultiVector<ST,LO,GO,NT>
51// thyra_X - destination Thyra::MultiVectorBase
52//
53void 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();
57
58 std::vector<RCP<Tpetra::MultiVector<ST, LO, GO, NT> > > subX;
59
60 // allocate vectors to copy into
61 Blocking::buildSubVectors(blockMaps_, subX, count);
62
63 // copy source vector to X vector
64 Blocking::one2many(subX, X, blockImport_);
65
66 // convert subX to an array of multi vectors
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);
74
75 fillDefaultSpmdMultiVector(vec, subX[i]);
76 }
77}
78
79// Virtual function defined in MappingStrategy. This copies
80// an Tpetra::MultiVector<ST,LO,GO,NT> into a Thyra::MultiVectorBase with
81// blocking handled by the strides defined in the constructor.
82//
83// arguments:
84// thyra_Y - source Thyra::MultiVectorBase
85// Y - destination Tpetra::MultiVector<ST,LO,GO,NT>
86//
87void 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);
93
94 // convert thyra product vector to subY
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());
100 }
101
102 // endow the subVectors with required information about the maps
103 // Blocking::associateSubVectors(blockMaps_,subY);
104
105 // copy solution vectors to Y vector
106 Blocking::many2one(Y, subY, blockExport_);
107}
108
109// this is the core routine that builds the maps
110// and importers/exporters neccessary for all the
111// transfers. Currently it simply calls out to the
112// interlaced tpetra functions. (Comment: this
113// routine should probably be private or protected
114// ... it is basically the meat of the constructor)
115//
116// arguments:
117// vars - Vector describing the blocking of variables
118// baseMap - basic map to use in the transfers
119// comm - Teuchos::RCP<Teuchos::Comm<int> > object
120//
121void 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) {
124 // build block for each vector
125 for (std::size_t i = 0; i < vars.size(); i++) {
126 // build maps and exporters/importers
127 Blocking::MapPair mapPair = Blocking::buildSubMap(vars[i], comm);
128 Blocking::ImExPair iePair = Blocking::buildExportImport(*baseMap, mapPair);
129
130 blockMaps_.push_back(mapPair);
131 blockImport_.push_back(iePair.first);
132 blockExport_.push_back(iePair.second);
133 }
134}
135
136// Builds a blocked Thyra operator that uses the strided
137// mapping strategy to define sub blocks.
138//
139// arguments:
140// mat - Tpetra::CrsMatrix<ST,LO,GO,NT> with FillComplete called, this
141// matrix is assumed to be square, with the same
142// range and domain maps
143// returns: Blocked Thyra linear operator with sub blocks
144// defined by this mapping strategy
145//
146const Teuchos::RCP<Thyra::BlockedLinearOpBase<ST> >
147TpetraBlockedMappingStrategy::buildBlockedThyraOp(
148 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
149 const std::string& label) const {
150 int dim = blockMaps_.size();
151
152 plocal2ContigGIDs.resize(dim);
153 for (int j = 0; j < dim; j++) {
154 plocal2ContigGIDs[j] = Blocking::getSubBlockColumnGIDs(*crsContent, blockMaps_[j]);
155 }
156
157 RCP<Thyra::DefaultBlockedLinearOp<ST> > A = Thyra::defaultBlockedLinearOp<ST>();
158
159 A->beginBlockFill(dim, dim);
160 for (int i = 0; i < dim; i++) {
161 for (int j = 0; j < dim; j++) {
162 // label block correctly
163 std::stringstream ss;
164 ss << label << "_" << i << "," << j;
165
166 // build the blocks and place it the right location
167 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > blk =
168 Blocking::buildSubBlock(i, j, crsContent, blockMaps_, plocal2ContigGIDs[j]);
169 A->setNonconstBlock(i, j,
170 Thyra::tpetraLinearOp<ST, LO, GO, NT>(
171 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(blk->getRangeMap()),
172 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(blk->getDomainMap()), blk));
173 }
174 } // end for i
175 A->endBlockFill();
176
177 return A;
178}
179
180// Rebuilds a blocked Thyra operator that uses the strided
181// mapping strategy to define sub blocks.
182//
183// arguments:
184// crsContent - Tpetra::CrsMatrix<ST,LO,GO,NT> with FillComplete called, this
185// matrix is assumed to be square, with the same
186// range and domain maps
187// A - Destination block linear op composed of blocks of
188// Tpetra::CrsMatrix<ST,LO,GO,NT> at all relevant locations
189//
190void TpetraBlockedMappingStrategy::rebuildBlockedThyraOp(
191 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
192 const RCP<Thyra::BlockedLinearOpBase<ST> >& A) const {
193 int dim = blockMaps_.size();
194 for (int i = 0; i < dim; i++) {
195 for (int j = 0; j < dim; j++) {
196 // get Tpetra version of desired block
197 RCP<Thyra::LinearOpBase<ST> > Aij = A->getNonconstBlock(i, j);
198 RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT> > tAij =
199 rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT> >(Aij, true);
200 RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > eAij =
201 rcp_dynamic_cast<Tpetra::CrsMatrix<ST, LO, GO, NT> >(tAij->getTpetraOperator(), true);
202
203 // rebuild the blocks and place it the right location
204 Blocking::rebuildSubBlock(i, j, crsContent, blockMaps_, *eAij, plocal2ContigGIDs[j]);
205 }
206 } // end for i
207}
208
209} // namespace TpetraHelpers
210} // end namespace Teko