00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047 #include "Teko_EpetraThyraConverter.hpp"
00048
00049
00050 #include "Teuchos_Array.hpp"
00051 #include "Teuchos_ArrayRCP.hpp"
00052
00053
00054 #include "Thyra_DefaultProductVectorSpace.hpp"
00055 #include "Thyra_DefaultProductMultiVector.hpp"
00056 #include "Thyra_SpmdMultiVectorBase.hpp"
00057 #include "Thyra_SpmdVectorSpaceBase.hpp"
00058 #include "Thyra_MultiVectorStdOps.hpp"
00059
00060 #include <iostream>
00061 #include <vector>
00062
00063 using Teuchos::RCP;
00064 using Teuchos::Ptr;
00065 using Teuchos::rcp;
00066 using Teuchos::rcpFromRef;
00067 using Teuchos::rcp_dynamic_cast;
00068 using Teuchos::ptr_dynamic_cast;
00069 using Teuchos::null;
00070
00071 namespace Teko {
00072 namespace Epetra {
00073
00074
00075
00076
00077 void blockEpetraToThyra(int numVectors,const double * epetraData,int leadingDim,const Teuchos::Ptr<Thyra::MultiVectorBase<double> > & mv,int & localDim)
00078 {
00079 localDim = 0;
00080
00081
00082 const Ptr<Thyra::ProductMultiVectorBase<double> > prodMV
00083 = ptr_dynamic_cast<Thyra::ProductMultiVectorBase<double> > (mv);
00084 if(prodMV==Teuchos::null) {
00085
00086 const Ptr<Thyra::SpmdMultiVectorBase<double> > spmdX = ptr_dynamic_cast<Thyra::SpmdMultiVectorBase<double> >(mv,true);
00087 const RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS = spmdX->spmdSpace();
00088
00089 int localSubDim = spmdVS->localSubDim();
00090
00091 Thyra::Ordinal thyraLeadingDim=0;
00092
00093
00094
00095 Teuchos::ArrayRCP<double> thyraData_arcp;
00096 Teuchos::ArrayView<double> thyraData;
00097 spmdX->getNonconstLocalData(Teuchos::outArg(thyraData_arcp),Teuchos::outArg(thyraLeadingDim));
00098 thyraData = thyraData_arcp();
00099
00100 for(int i=0;i<localSubDim;i++) {
00101
00102 for(int v=0;v<numVectors;v++)
00103 thyraData[i+thyraLeadingDim*v] = epetraData[i+leadingDim*v];
00104 }
00105
00106
00107
00108
00109 localDim = localSubDim;
00110
00111 return;
00112 }
00113
00114
00115 const double * localData = epetraData;
00116
00117
00118 for(int blkIndex=0;blkIndex<prodMV->productSpace()->numBlocks();blkIndex++) {
00119 int subDim = 0;
00120 const RCP<Thyra::MultiVectorBase<double> > blockVec = prodMV->getNonconstMultiVectorBlock(blkIndex);
00121
00122
00123 blockEpetraToThyra(numVectors, localData,leadingDim,blockVec.ptr(),subDim);
00124
00125
00126 localData += subDim;
00127
00128
00129 localDim += subDim;
00130 }
00131 }
00132
00133
00134
00135
00136 void blockEpetraToThyra(const Epetra_MultiVector & epetraX,const Teuchos::Ptr<Thyra::MultiVectorBase<double> > & thyraX)
00137 {
00138 TEUCHOS_ASSERT(thyraX->range()->dim()==epetraX.GlobalLength());
00139
00140
00141 int leadingDim=0,numVectors=0,localDim=0;
00142 double * epetraData=0;
00143 epetraX.ExtractView(&epetraData,&leadingDim);
00144
00145 numVectors = epetraX.NumVectors();
00146
00147 blockEpetraToThyra(numVectors,epetraData,leadingDim,thyraX.ptr(),localDim);
00148
00149 TEUCHOS_ASSERT(localDim==epetraX.MyLength());
00150 }
00151
00152 void blockThyraToEpetra(int numVectors,double * epetraData,int leadingDim,const Teuchos::RCP<const Thyra::MultiVectorBase<double> > & tX,int & localDim)
00153 {
00154 localDim = 0;
00155
00156
00157 const RCP<const Thyra::ProductMultiVectorBase<double> > prodX
00158 = rcp_dynamic_cast<const Thyra::ProductMultiVectorBase<double> > (tX);
00159 if(prodX==Teuchos::null) {
00160
00161
00162
00163 RCP<const Thyra::SpmdMultiVectorBase<double> > spmdX = rcp_dynamic_cast<const Thyra::SpmdMultiVectorBase<double> >(tX,true);
00164 RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS = spmdX->spmdSpace();
00165
00166 int localSubDim = spmdVS->localSubDim();
00167
00168 Thyra::Ordinal thyraLeadingDim=0;
00169
00170
00171
00172 Teuchos::ArrayView<const double> thyraData;
00173 Teuchos::ArrayRCP<const double> thyraData_arcp;
00174 spmdX->getLocalData(Teuchos::outArg(thyraData_arcp),Teuchos::outArg(thyraLeadingDim));
00175 thyraData = thyraData_arcp();
00176
00177 for(int i=0;i<localSubDim;i++) {
00178
00179 for(int v=0;v<numVectors;v++)
00180 epetraData[i+leadingDim*v] = thyraData[i+thyraLeadingDim*v];
00181 }
00182
00183
00184 localDim = localSubDim;
00185
00186 return;
00187 }
00188
00189 const RCP<const Thyra::ProductVectorSpaceBase<double> > prodVS = prodX->productSpace();
00190
00191
00192 double * localData = epetraData;
00193
00194
00195 for(int blkIndex=0;blkIndex<prodVS->numBlocks();blkIndex++) {
00196 int subDim = 0;
00197
00198
00199 blockThyraToEpetra(numVectors, localData,leadingDim,prodX->getMultiVectorBlock(blkIndex),subDim);
00200
00201
00202 localData += subDim;
00203
00204
00205 localDim += subDim;
00206 }
00207
00208 return;
00209 }
00210
00211
00212
00213
00214
00215 void blockThyraToEpetra(const Teuchos::RCP<const Thyra::MultiVectorBase<double> > & thyraX,Epetra_MultiVector & epetraX)
00216 {
00217
00218 int numVectors = thyraX->domain()->dim();
00219
00220
00221 TEUCHOS_ASSERT(numVectors==epetraX.NumVectors());
00222 TEUCHOS_ASSERT(thyraX->range()->dim()==epetraX.GlobalLength());
00223
00224
00225 int leadingDim=0,localDim=0;
00226 double * epetraData=0;
00227 epetraX.ExtractView(&epetraData,&leadingDim);
00228
00229
00230 blockThyraToEpetra(numVectors,epetraData,leadingDim,thyraX,localDim);
00231
00232
00233 TEUCHOS_ASSERT(localDim==epetraX.Map().NumMyElements());
00234 }
00235
00236 void thyraVSToEpetraMap(std::vector<int> & myIndicies, int blockOffset, const Thyra::VectorSpaceBase<double> & vs, int & localDim)
00237 {
00238
00239 localDim = 0;
00240
00241 const RCP<const Thyra::ProductVectorSpaceBase<double> > prodVS
00242 = rcp_dynamic_cast<const Thyra::ProductVectorSpaceBase<double> >(rcpFromRef(vs));
00243
00244
00245 if(prodVS==Teuchos::null) {
00246
00247
00248
00249 const RCP<const Thyra::SpmdVectorSpaceBase<double> > spmdVS
00250 = rcp_dynamic_cast<const Thyra::SpmdVectorSpaceBase<double> >(rcpFromRef(vs));
00251 TEUCHOS_TEST_FOR_EXCEPTION(spmdVS==Teuchos::null,std::runtime_error,
00252 "thyraVSToEpetraMap requires all subblocks to be SPMD");
00253
00254
00255 int localOffset = spmdVS->localOffset();
00256 int localSubDim = spmdVS->localSubDim();
00257
00258
00259 for(int i=0;i<localSubDim;i++)
00260 myIndicies.push_back(blockOffset+localOffset+i);
00261
00262 localDim += localSubDim;
00263
00264 return;
00265 }
00266
00267
00268 for(int blkIndex=0;blkIndex<prodVS->numBlocks();blkIndex++) {
00269 int subDim = 0;
00270
00271
00272 thyraVSToEpetraMap(myIndicies, blockOffset,*prodVS->getBlock(blkIndex),subDim);
00273
00274 blockOffset += prodVS->getBlock(blkIndex)->dim();
00275
00276
00277 localDim += subDim;
00278 }
00279 }
00280
00281
00282 const RCP<Epetra_Map> thyraVSToEpetraMap(const Thyra::VectorSpaceBase<double> & vs,const RCP<const Epetra_Comm> & comm)
00283 {
00284 int localDim = 0;
00285 std::vector<int> myGIDs;
00286
00287
00288 thyraVSToEpetraMap(myGIDs,0,vs,localDim);
00289
00290 TEUCHOS_ASSERT(myGIDs.size()==(unsigned int) localDim);
00291
00292
00293 return rcp(new Epetra_Map(vs.dim(), myGIDs.size(), &(myGIDs[0]), 0, *comm));
00294 }
00295
00296 }
00297 }