Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_TpetraStridedMappingStrategy.cpp
1 /*
2 // @HEADER
3 //
4 // ***********************************************************************
5 //
6 // Teko: A package for block and physics based preconditioning
7 // Copyright 2010 Sandia Corporation
8 //
9 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
10 // the U.S. Government retains certain rights in this software.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
40 //
41 // ***********************************************************************
42 //
43 // @HEADER
44 
45 */
46 
47 #include "Teko_TpetraStridedMappingStrategy.hpp"
48 #include "Teko_InterlacedTpetra.hpp"
49 #include "Teko_TpetraHelpers.hpp"
50 
51 #include "Thyra_TpetraThyraWrappers.hpp"
52 #include "Thyra_TpetraLinearOp.hpp"
53 #include "Thyra_DefaultProductMultiVector.hpp"
54 #include "Thyra_DefaultProductVectorSpace.hpp"
55 #include "Thyra_DefaultSpmdMultiVector.hpp"
56 #include "Thyra_DefaultBlockedLinearOp.hpp"
57 
58 using Teuchos::RCP;
59 using Teuchos::rcp;
60 using Teuchos::rcp_dynamic_cast;
61 
62 namespace Teko {
63 namespace TpetraHelpers {
64 
65 // Creates a strided mapping strategy. This class is useful
66 // for breaking up nodally ordered matrices (i.e. the unknowns
67 // in a FEM problem are ordered [u0,v0,p0,u1,v1,p1,...]). Current
68 // implimentation only supports a fixed number of variables
69 //
70 // arguments:
71 // vars - Number of different variables
72 // map - original Tpetra::Map<LO,GO,NT> to be broken up
73 // comm - Teuchos::Comm<int> object related to the map
74 //
75 TpetraStridedMappingStrategy::TpetraStridedMappingStrategy(
76  const std::vector<int>& vars, const RCP<const Tpetra::Map<LO, GO, NT> >& map,
77  const Teuchos::Comm<int>& comm) {
78  rangeMap_ = map;
79  domainMap_ = map;
80  buildBlockTransferData(vars, rangeMap_, comm);
81 }
82 
83 // Virtual function defined in MappingStrategy. This copies
84 // an Tpetra::MultiVector<ST,LO,GO,NT> into a Thyra::MultiVectorBase with
85 // blocking handled by the strides defined in the constructor.
86 //
87 // arguments:
88 // X - source Tpetra::MultiVector<ST,LO,GO,NT>
89 // thyra_X - destination Thyra::MultiVectorBase
90 //
91 void TpetraStridedMappingStrategy::copyTpetraIntoThyra(
92  const Tpetra::MultiVector<ST, LO, GO, NT>& X,
93  const Teuchos::Ptr<Thyra::MultiVectorBase<ST> >& thyra_X) const {
94  int count = X.getNumVectors();
95 
96  std::vector<RCP<Tpetra::MultiVector<ST, LO, GO, NT> > > subX;
97 
98  // allocate vectors to copy into
99  Strided::buildSubVectors(blockMaps_, subX, count);
100 
101  // copy source vector to X vector
102  Strided::one2many(subX, X, blockImport_);
103 
104  // convert subX to an array of multi vectors
105  Teuchos::Ptr<Thyra::ProductMultiVectorBase<ST> > prod_X =
106  Teuchos::ptr_dynamic_cast<Thyra::ProductMultiVectorBase<ST> >(thyra_X);
107  for (unsigned int i = 0; i < blockMaps_.size(); i++) {
108  RCP<Thyra::TpetraMultiVector<ST, LO, GO, NT> > vec =
109  rcp_dynamic_cast<Thyra::TpetraMultiVector<ST, LO, GO, NT> >(
110  prod_X->getNonconstMultiVectorBlock(i), true);
111 
112  fillDefaultSpmdMultiVector(vec, subX[i]);
113  }
114 }
115 
116 // Virtual function defined in MappingStrategy. This copies
117 // an Tpetra::MultiVector<ST,LO,GO,NT> into a Thyra::MultiVectorBase with
118 // blocking handled by the strides defined in the constructor.
119 //
120 // arguments:
121 // thyra_Y - source Thyra::MultiVectorBase
122 // Y - destination Tpetra::MultiVector<ST,LO,GO,NT>
123 //
124 void TpetraStridedMappingStrategy::copyThyraIntoTpetra(
125  const RCP<const Thyra::MultiVectorBase<ST> >& thyra_Y,
126  Tpetra::MultiVector<ST, LO, GO, NT>& Y) const {
127  std::vector<RCP<const Tpetra::MultiVector<ST, LO, GO, NT> > > subY;
128  RCP<const Thyra::DefaultProductMultiVector<ST> > prod_Y =
129  rcp_dynamic_cast<const Thyra::DefaultProductMultiVector<ST> >(thyra_Y);
130 
131  // convert thyra product vector to subY
132  for (unsigned int i = 0; i < blockMaps_.size(); i++) {
133  RCP<const Thyra::TpetraMultiVector<ST, LO, GO, NT> > tmv =
134  rcp_dynamic_cast<const Thyra::TpetraMultiVector<ST, LO, GO, NT> >(
135  prod_Y->getMultiVectorBlock(i), true);
136  subY.push_back(tmv->getConstTpetraMultiVector());
137  }
138 
139  // endow the subVectors with required information about the maps
140  Strided::associateSubVectors(blockMaps_, subY);
141 
142  // copy solution vectors to Y vector
143  Strided::many2one(Y, subY, blockExport_);
144 }
145 
146 // this is the core routine that builds the maps
147 // and importers/exporters neccessary for all the
148 // transfers. Currently it simply calls out to the
149 // interlaced tpetra functions. (Comment: this
150 // routine should probably be private or protected
151 // ... it is basically the meat of the constructor)
152 //
153 // arguments:
154 // vars - Vector describing the blocking of variables
155 // baseMap - basic map to use in the transfers
156 // comm - Teuchos::Comm<int> object
157 //
158 void TpetraStridedMappingStrategy::buildBlockTransferData(
159  const std::vector<int>& vars, const Teuchos::RCP<const Tpetra::Map<LO, GO, NT> >& baseMap,
160  const Teuchos::Comm<int>& comm) {
161  // build maps and exporters/importers
162  Strided::buildSubMaps(*baseMap, vars, comm, blockMaps_);
163  Strided::buildExportImport(*baseMap, blockMaps_, blockExport_, blockImport_);
164 }
165 
166 // Builds a blocked Thyra operator that uses the strided
167 // mapping strategy to define sub blocks.
168 //
169 // arguments:
170 // mat - Tpetra::CrsMatrix<ST,LO,GO,NT> with FillComplete called, this
171 // matrix is assumed to be square, with the same
172 // range and domain maps
173 // returns: Blocked Thyra linear operator with sub blocks
174 // defined by this mapping strategy
175 //
176 const Teuchos::RCP<Thyra::BlockedLinearOpBase<ST> >
177 TpetraStridedMappingStrategy::buildBlockedThyraOp(
178  const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
179  const std::string& label) const {
180  int dim = blockMaps_.size();
181 
182  RCP<Thyra::DefaultBlockedLinearOp<ST> > A = Thyra::defaultBlockedLinearOp<ST>();
183 
184  A->beginBlockFill(dim, dim);
185  for (int i = 0; i < dim; i++) {
186  for (int j = 0; j < dim; j++) {
187  // label block correctly
188  std::stringstream ss;
189  ss << label << "_" << i << "," << j;
190 
191  // build the blocks and place it the right location
192  RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > Aij =
193  Strided::buildSubBlock(i, j, crsContent, blockMaps_);
194  A->setNonconstBlock(i, j,
195  Thyra::tpetraLinearOp<ST, LO, GO, NT>(
196  Thyra::tpetraVectorSpace<ST, LO, GO, NT>(Aij->getRangeMap()),
197  Thyra::tpetraVectorSpace<ST, LO, GO, NT>(Aij->getDomainMap()), Aij));
198  }
199  } // end for i
200  A->endBlockFill();
201 
202  return A;
203 }
204 
205 // Rebuilds a blocked Thyra operator that uses the strided
206 // mapping strategy to define sub blocks.
207 //
208 // arguments:
209 // crsContent - Tpetra::CrsMatrix<ST,LO,GO,NT> with FillComplete called, this
210 // matrix is assumed to be square, with the same
211 // range and domain maps
212 // A - Destination block linear op composed of blocks of
213 // Tpetra::CrsMatrix<ST,LO,GO,NT> at all relevant locations
214 //
215 void TpetraStridedMappingStrategy::rebuildBlockedThyraOp(
216  const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> >& crsContent,
217  const RCP<Thyra::BlockedLinearOpBase<ST> >& A) const {
218  int dim = blockMaps_.size();
219 
220  for (int i = 0; i < dim; i++) {
221  for (int j = 0; j < dim; j++) {
222  // get Tpetra version of desired block
223  RCP<Thyra::LinearOpBase<ST> > Aij = A->getNonconstBlock(i, j);
224  RCP<Thyra::TpetraLinearOp<ST, LO, GO, NT> > tAij =
225  rcp_dynamic_cast<Thyra::TpetraLinearOp<ST, LO, GO, NT> >(Aij, true);
226  RCP<Tpetra::CrsMatrix<ST, LO, GO, NT> > eAij =
227  rcp_dynamic_cast<Tpetra::CrsMatrix<ST, LO, GO, NT> >(tAij->getTpetraOperator(), true);
228 
229  // rebuild the blocks and place it the right location
230  Strided::rebuildSubBlock(i, j, crsContent, blockMaps_, *eAij);
231  }
232  } // end for i
233 }
234 
235 } // end namespace TpetraHelpers
236 } // end namespace Teko