47 #include "Teko_EpetraThyraConverter.hpp"
50 #include "Teuchos_Array.hpp"
51 #include "Teuchos_ArrayRCP.hpp"
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"
66 using Teuchos::rcpFromRef;
67 using Teuchos::rcp_dynamic_cast;
68 using Teuchos::ptr_dynamic_cast;
77 void blockEpetraToThyra(
int numVectors,
const double * epetraData,
int leadingDim,
const Teuchos::Ptr<Thyra::MultiVectorBase<double> > & mv,
int & localDim)
82 const Ptr<Thyra::ProductMultiVectorBase<double> > prodMV
83 = ptr_dynamic_cast<Thyra::ProductMultiVectorBase<double> > (mv);
84 if(prodMV==Teuchos::null) {
86 const Ptr<Thyra::SpmdMultiVectorBase<double> > spmdX = ptr_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(mv,
true);
87 const RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS = spmdX->spmdSpace();
89 int localSubDim = spmdVS->localSubDim();
91 Thyra::Ordinal thyraLeadingDim=0;
95 Teuchos::ArrayRCP<double> thyraData_arcp;
96 Teuchos::ArrayView<double> thyraData;
97 spmdX->getNonconstLocalData(Teuchos::outArg(thyraData_arcp),Teuchos::outArg(thyraLeadingDim));
98 thyraData = thyraData_arcp();
100 for(
int i=0;i<localSubDim;i++) {
102 for(
int v=0;v<numVectors;v++)
103 thyraData[i+thyraLeadingDim*v] = epetraData[i+leadingDim*v];
109 localDim = localSubDim;
115 const double * localData = epetraData;
118 for(
int blkIndex=0;blkIndex<prodMV->productSpace()->numBlocks();blkIndex++) {
120 const RCP<Thyra::MultiVectorBase<double> > blockVec = prodMV->getNonconstMultiVectorBlock(blkIndex);
123 blockEpetraToThyra(numVectors, localData,leadingDim,blockVec.ptr(),subDim);
136 void blockEpetraToThyra(
const Epetra_MultiVector & epetraX,
const Teuchos::Ptr<Thyra::MultiVectorBase<double> > & thyraX)
138 TEUCHOS_ASSERT(thyraX->range()->dim()==epetraX.GlobalLength());
141 int leadingDim=0,numVectors=0,localDim=0;
142 double * epetraData=0;
143 epetraX.ExtractView(&epetraData,&leadingDim);
145 numVectors = epetraX.NumVectors();
147 blockEpetraToThyra(numVectors,epetraData,leadingDim,thyraX.ptr(),localDim);
149 TEUCHOS_ASSERT(localDim==epetraX.MyLength());
152 void blockThyraToEpetra(
int numVectors,
double * epetraData,
int leadingDim,
const Teuchos::RCP<
const Thyra::MultiVectorBase<double> > & tX,
int & localDim)
157 const RCP<const Thyra::ProductMultiVectorBase<double> > prodX
158 = rcp_dynamic_cast<
const Thyra::ProductMultiVectorBase<double> > (tX);
159 if(prodX==Teuchos::null) {
163 RCP<const Thyra::SpmdMultiVectorBase<double> > spmdX = rcp_dynamic_cast<
const Thyra::SpmdMultiVectorBase<double> >(tX,
true);
164 RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS = spmdX->spmdSpace();
166 int localSubDim = spmdVS->localSubDim();
168 Thyra::Ordinal thyraLeadingDim=0;
172 Teuchos::ArrayView<const double> thyraData;
173 Teuchos::ArrayRCP<const double> thyraData_arcp;
174 spmdX->getLocalData(Teuchos::outArg(thyraData_arcp),Teuchos::outArg(thyraLeadingDim));
175 thyraData = thyraData_arcp();
177 for(
int i=0;i<localSubDim;i++) {
179 for(
int v=0;v<numVectors;v++)
180 epetraData[i+leadingDim*v] = thyraData[i+thyraLeadingDim*v];
184 localDim = localSubDim;
189 const RCP<const Thyra::ProductVectorSpaceBase<double> > prodVS = prodX->productSpace();
192 double * localData = epetraData;
195 for(
int blkIndex=0;blkIndex<prodVS->numBlocks();blkIndex++) {
199 blockThyraToEpetra(numVectors, localData,leadingDim,prodX->getMultiVectorBlock(blkIndex),subDim);
215 void blockThyraToEpetra(
const Teuchos::RCP<
const Thyra::MultiVectorBase<double> > & thyraX,Epetra_MultiVector & epetraX)
218 int numVectors = thyraX->domain()->dim();
221 TEUCHOS_ASSERT(numVectors==epetraX.NumVectors());
222 TEUCHOS_ASSERT(thyraX->range()->dim()==epetraX.GlobalLength());
225 int leadingDim=0,localDim=0;
226 double * epetraData=0;
227 epetraX.ExtractView(&epetraData,&leadingDim);
230 blockThyraToEpetra(numVectors,epetraData,leadingDim,thyraX,localDim);
233 TEUCHOS_ASSERT(localDim==epetraX.Map().NumMyElements());
236 void thyraVSToEpetraMap(std::vector<int> & myIndicies,
int blockOffset,
const Thyra::VectorSpaceBase<double> & vs,
int & localDim)
241 const RCP<const Thyra::ProductVectorSpaceBase<double> > prodVS
242 = rcp_dynamic_cast<
const Thyra::ProductVectorSpaceBase<double> >(rcpFromRef(vs));
245 if(prodVS==Teuchos::null) {
249 const RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS
250 = rcp_dynamic_cast<
const Thyra::SpmdVectorSpaceBase<double> >(rcpFromRef(vs));
251 TEUCHOS_TEST_FOR_EXCEPTION(spmdVS==Teuchos::null,std::runtime_error,
252 "thyraVSToEpetraMap requires all subblocks to be SPMD");
255 int localOffset = spmdVS->localOffset();
256 int localSubDim = spmdVS->localSubDim();
259 for(
int i=0;i<localSubDim;i++)
260 myIndicies.push_back(blockOffset+localOffset+i);
262 localDim += localSubDim;
268 for(
int blkIndex=0;blkIndex<prodVS->numBlocks();blkIndex++) {
272 thyraVSToEpetraMap(myIndicies, blockOffset,*prodVS->getBlock(blkIndex),subDim);
274 blockOffset += prodVS->getBlock(blkIndex)->dim();
282 const RCP<Epetra_Map> thyraVSToEpetraMap(
const Thyra::VectorSpaceBase<double> & vs,
const RCP<const Epetra_Comm> & comm)
285 std::vector<int> myGIDs;
288 thyraVSToEpetraMap(myGIDs,0,vs,localDim);
290 TEUCHOS_ASSERT(myGIDs.size()==(
unsigned int) localDim);
293 return rcp(
new Epetra_Map(vs.dim(), myGIDs.size(), &(myGIDs[0]), 0, *comm));