Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_TpetraThyraConverter.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_TpetraThyraConverter.hpp"
48 #include "Tpetra_Core.hpp"
49 
50 // Teuchos includes
51 #include "Teuchos_Array.hpp"
52 #include "Teuchos_ArrayRCP.hpp"
53 
54 // Thyra includes
55 #include "Thyra_DefaultProductVectorSpace.hpp"
56 #include "Thyra_DefaultProductMultiVector.hpp"
57 #include "Thyra_SpmdMultiVectorBase.hpp"
58 #include "Thyra_SpmdVectorSpaceBase.hpp"
59 #include "Thyra_MultiVectorStdOps.hpp"
60 
61 #include <iostream>
62 #include <vector>
63 
64 using Teuchos::null;
65 using Teuchos::Ptr;
66 using Teuchos::ptr_dynamic_cast;
67 using Teuchos::RCP;
68 using Teuchos::rcp;
69 using Teuchos::rcp_dynamic_cast;
70 using Teuchos::rcpFromRef;
71 
72 namespace Teko {
73 namespace TpetraHelpers {
74 
75 // const Teuchos::RCP<const Thyra::MultiVectorBase<double> >
76 // blockTpetraToThyra(int numVectors,const double * tpetraData,int leadingDim,const
77 // Teuchos::RCP<const Thyra::VectorSpaceBase<double> > & vs,int & localDim)
78 
79 void blockTpetraToThyra(int numVectors, Teuchos::ArrayRCP<const ST> tpetraData, int leadingDim,
80  const Teuchos::Ptr<Thyra::MultiVectorBase<ST> >& mv, int& localDim) {
81  localDim = 0;
82 
83  // check the base case
84  const Ptr<Thyra::ProductMultiVectorBase<ST> > prodMV =
85  ptr_dynamic_cast<Thyra::ProductMultiVectorBase<ST> >(mv);
86  if (prodMV == Teuchos::null) {
87  // VS object must be a SpmdMultiVector object
88  const Ptr<Thyra::SpmdMultiVectorBase<ST> > spmdX =
89  ptr_dynamic_cast<Thyra::SpmdMultiVectorBase<ST> >(mv, true);
90  const RCP<const Thyra::SpmdVectorSpaceBase<ST> > spmdVS = spmdX->spmdSpace();
91 
92  int localSubDim = spmdVS->localSubDim();
93 
94  Thyra::Ordinal thyraLeadingDim = 0;
95 
96  Teuchos::ArrayRCP<ST> thyraData_arcp;
97  Teuchos::ArrayView<ST> thyraData;
98  spmdX->getNonconstLocalData(Teuchos::outArg(thyraData_arcp), Teuchos::outArg(thyraLeadingDim));
99  thyraData = thyraData_arcp(); // build array view
100 
101  for (int i = 0; i < localSubDim; i++) {
102  // copy each vector
103  for (int v = 0; v < numVectors; v++)
104  thyraData[i + thyraLeadingDim * v] = tpetraData[i + leadingDim * v];
105  }
106 
107  // set the local dimension
108  localDim = localSubDim;
109 
110  return;
111  }
112 
113  // this keeps track of current location in the tpetraData vector
114  Teuchos::ArrayRCP<const ST> localData = tpetraData;
115 
116  // loop over all the blocks in the vector space
117  for (int blkIndex = 0; blkIndex < prodMV->productSpace()->numBlocks(); blkIndex++) {
118  int subDim = 0;
119  const RCP<Thyra::MultiVectorBase<ST> > blockVec = prodMV->getNonconstMultiVectorBlock(blkIndex);
120 
121  // perform the recusive copy
122  blockTpetraToThyra(numVectors, localData, leadingDim, blockVec.ptr(), subDim);
123 
124  // shift to the next block
125  localData += subDim;
126 
127  // account for the size of this subblock
128  localDim += subDim;
129  }
130 }
131 
132 // Convert a Tpetra_MultiVector with assumed block structure dictated by the
133 // vector space into a Thyra::MultiVectorBase object.
134 // const Teuchos::RCP<const Thyra::MultiVectorBase<double> > blockTpetraToThyra(const
135 // Tpetra_MultiVector & e,const Teuchos::RCP<const Thyra::VectorSpaceBase<double> > & vs)
136 void blockTpetraToThyra(const Tpetra::MultiVector<ST, LO, GO, NT>& tpetraX,
137  const Teuchos::Ptr<Thyra::MultiVectorBase<ST> >& thyraX) {
138  TEUCHOS_ASSERT((Tpetra::global_size_t)thyraX->range()->dim() == tpetraX.getGlobalLength());
139 
140  // extract local information from the Tpetra_MultiVector
141  LO leadingDim = 0, localDim = 0;
142  leadingDim = tpetraX.getStride();
143  Teuchos::ArrayRCP<const ST> tpetraData = tpetraX.get1dView();
144 
145  int numVectors = tpetraX.getNumVectors();
146 
147  blockTpetraToThyra(numVectors, tpetraData, leadingDim, thyraX.ptr(), localDim);
148 
149  TEUCHOS_ASSERT((size_t)localDim == tpetraX.getLocalLength());
150 }
151 
152 void blockThyraToTpetra(LO numVectors, Teuchos::ArrayRCP<ST> tpetraData, LO leadingDim,
153  const Teuchos::RCP<const Thyra::MultiVectorBase<ST> >& tX, LO& localDim) {
154  localDim = 0;
155 
156  // check the base case
157  const RCP<const Thyra::ProductMultiVectorBase<ST> > prodX =
158  rcp_dynamic_cast<const Thyra::ProductMultiVectorBase<ST> >(tX);
159  if (prodX == Teuchos::null) {
160  // the base case
161 
162  // VS object must be a SpmdMultiVector object
163  RCP<const Thyra::SpmdMultiVectorBase<ST> > spmdX =
164  rcp_dynamic_cast<const Thyra::SpmdMultiVectorBase<ST> >(tX, true);
165  RCP<const Thyra::SpmdVectorSpaceBase<ST> > spmdVS = spmdX->spmdSpace();
166 
167  Thyra::Ordinal thyraLeadingDim = 0;
168  Teuchos::ArrayView<const ST> thyraData;
169  Teuchos::ArrayRCP<const ST> thyraData_arcp;
170  spmdX->getLocalData(Teuchos::outArg(thyraData_arcp), Teuchos::outArg(thyraLeadingDim));
171  thyraData = thyraData_arcp(); // grab the array view
172 
173  LO localSubDim = spmdVS->localSubDim();
174  for (LO i = 0; i < localSubDim; i++) {
175  // copy each vector
176  for (LO v = 0; v < numVectors; v++) {
177  tpetraData[i + leadingDim * v] = thyraData[i + thyraLeadingDim * v];
178  }
179  }
180 
181  // set the local dimension
182  localDim = localSubDim;
183 
184  return;
185  }
186 
187  const RCP<const Thyra::ProductVectorSpaceBase<ST> > prodVS = prodX->productSpace();
188 
189  // this keeps track of current location in the tpetraData vector
190  Teuchos::ArrayRCP<ST> localData = tpetraData;
191 
192  // loop over all the blocks in the vector space
193  for (int blkIndex = 0; blkIndex < prodVS->numBlocks(); blkIndex++) {
194  int subDim = 0;
195 
196  // construct the block vector
197  blockThyraToTpetra(numVectors, localData, leadingDim, prodX->getMultiVectorBlock(blkIndex),
198  subDim);
199 
200  // shift to the next block
201  localData += subDim;
202 
203  // account for the size of this subblock
204  localDim += subDim;
205  }
206 
207  return;
208 }
209 
210 // Convert a Thyra::MultiVectorBase object to a Tpetra_MultiVector object with
211 // the map defined by the Tpetra_Map.
212 // const Teuchos::RCP<const Tpetra_MultiVector>
213 // blockThyraToTpetra(const Teuchos::RCP<const Thyra::MultiVectorBase<double> > & tX,const RCP<const
214 // Tpetra_Map> & map)
215 void blockThyraToTpetra(const Teuchos::RCP<const Thyra::MultiVectorBase<ST> >& thyraX,
216  Tpetra::MultiVector<ST, LO, GO, NT>& tpetraX) {
217  // build an Tpetra_MultiVector object
218  LO numVectors = thyraX->domain()->dim();
219 
220  // make sure the number of vectors are the same
221  TEUCHOS_ASSERT((size_t)numVectors == tpetraX.getNumVectors());
222  TEUCHOS_ASSERT((Tpetra::global_size_t)thyraX->range()->dim() == tpetraX.getGlobalLength());
223 
224  // extract local information from the Tpetra_MultiVector
225  LO leadingDim = 0, localDim = 0;
226  leadingDim = tpetraX.getStride();
227  Teuchos::ArrayRCP<ST> tpetraData = tpetraX.get1dViewNonConst();
228 
229  // perform recursive copy
230  blockThyraToTpetra(numVectors, tpetraData, leadingDim, thyraX, localDim);
231 
232  // sanity check
233  TEUCHOS_ASSERT((size_t)localDim == tpetraX.getLocalLength());
234 }
235 
236 void thyraVSToTpetraMap(std::vector<GO>& myIndicies, int blockOffset,
237  const Thyra::VectorSpaceBase<ST>& vs, int& localDim) {
238  // zero out set local dimension
239  localDim = 0;
240 
241  const RCP<const Thyra::ProductVectorSpaceBase<ST> > prodVS =
242  rcp_dynamic_cast<const Thyra::ProductVectorSpaceBase<ST> >(rcpFromRef(vs));
243 
244  // is more recursion needed?
245  if (prodVS == Teuchos::null) {
246  // base case
247 
248  // try to cast to an SPMD capable vector space
249  const RCP<const Thyra::SpmdVectorSpaceBase<ST> > spmdVS =
250  rcp_dynamic_cast<const Thyra::SpmdVectorSpaceBase<ST> >(rcpFromRef(vs));
251  TEUCHOS_TEST_FOR_EXCEPTION(spmdVS == Teuchos::null, std::runtime_error,
252  "thyraVSToTpetraMap requires all subblocks to be SPMD");
253 
254  // get local data storage information
255  int localOffset = spmdVS->localOffset();
256  int localSubDim = spmdVS->localSubDim();
257 
258  // add indicies to matrix
259  for (int i = 0; i < localSubDim; i++) myIndicies.push_back(blockOffset + localOffset + i);
260 
261  localDim += localSubDim;
262 
263  return;
264  }
265 
266  // loop over all the blocks in the vector space
267  for (int blkIndex = 0; blkIndex < prodVS->numBlocks(); blkIndex++) {
268  int subDim = 0;
269 
270  // construct the block vector
271  thyraVSToTpetraMap(myIndicies, blockOffset, *prodVS->getBlock(blkIndex), subDim);
272 
273  blockOffset += prodVS->getBlock(blkIndex)->dim();
274 
275  // account for the size of this subblock
276  localDim += subDim;
277  }
278 }
279 
280 // From a Thyra vector space create a compatable Tpetra_Map
281 const RCP<Tpetra::Map<LO, GO, NT> > thyraVSToTpetraMap(
282  const Thyra::VectorSpaceBase<ST>& vs,
283  const RCP<const Teuchos::Comm<Thyra::Ordinal> >& /* comm */) {
284  int localDim = 0;
285  std::vector<GO> myGIDs;
286 
287  // call recursive routine that constructs the mapping
288  thyraVSToTpetraMap(myGIDs, 0, vs, localDim);
289 
290  TEUCHOS_ASSERT(myGIDs.size() == (size_t)localDim);
291 
292  // FIXME (mfh 12 Jul 2018) This ignores the input comm, so it can't
293  // be right.
294 
295  // create the map
296  return rcp(new Tpetra::Map<LO, GO, NT>(vs.dim(), Teuchos::ArrayView<const GO>(myGIDs), 0,
297  Tpetra::getDefaultComm()));
298 }
299 
300 } // namespace TpetraHelpers
301 } // end namespace Teko