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