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_BlockingEpetra.hpp"
00048
00049 #include "Epetra_IntVector.h"
00050 #include "Epetra_LocalMap.h"
00051
00052 using Teuchos::RCP;
00053 using Teuchos::rcp;
00054
00055 namespace Teko {
00056 namespace Epetra {
00057 namespace Blocking {
00058
00072 const MapPair buildSubMap(const std::vector< int > & gid, const Epetra_Comm &comm)
00073 {
00074 Teuchos::RCP<Epetra_Map> gidMap = rcp(new Epetra_Map(-1,gid.size(),&gid[0],0,comm));
00075 Teuchos::RCP<Epetra_Map> contigMap = rcp(new Epetra_Map(-1,gid.size(),0,comm));
00076
00077 return std::make_pair(gidMap,contigMap);
00078 }
00079
00089 const ImExPair buildExportImport(const Epetra_Map & baseMap,const MapPair & maps)
00090 {
00091 return std::make_pair(rcp(new Epetra_Import(*maps.first,baseMap)),
00092 rcp(new Epetra_Export(*maps.first,baseMap)));
00093 }
00094
00102 void buildSubVectors(const std::vector<MapPair> & maps,
00103 std::vector<RCP<Epetra_MultiVector> > & vectors,int count)
00104 {
00105 std::vector<MapPair>::const_iterator mapItr;
00106
00107
00108 for(mapItr=maps.begin();mapItr!=maps.end();mapItr++) {
00109
00110 RCP<Epetra_MultiVector> mv = rcp(new Epetra_MultiVector(*(*mapItr).second,count));
00111 vectors.push_back(mv);
00112 }
00113 }
00114
00121 void one2many(std::vector<RCP<Epetra_MultiVector> > & many, const Epetra_MultiVector & single,
00122 const std::vector<RCP<Epetra_Import> > & subImport)
00123 {
00124
00125 std::vector<RCP<Epetra_MultiVector> >::const_iterator vecItr;
00126 std::vector<RCP<Epetra_Import> >::const_iterator impItr;
00127
00128
00129 for(vecItr=many.begin(),impItr=subImport.begin();
00130 vecItr!=many.end();++vecItr,++impItr) {
00131
00132 RCP<Epetra_MultiVector> destVec = *vecItr;
00133
00134
00135 const Epetra_BlockMap & globalMap = (*impItr)->TargetMap();
00136
00137
00138 Epetra_MultiVector importVector(View,globalMap,destVec->Values(),destVec->Stride(),destVec->NumVectors());
00139
00140
00141 importVector.Import(single,**impItr,Insert);
00142 }
00143 }
00144
00154 void many2one(Epetra_MultiVector & one, const std::vector<RCP<const Epetra_MultiVector> > & many,
00155 const std::vector<RCP<Epetra_Export> > & subExport)
00156 {
00157
00158 std::vector<RCP<const Epetra_MultiVector> >::const_iterator vecItr;
00159 std::vector<RCP<Epetra_Export> >::const_iterator expItr;
00160
00161
00162 for(vecItr=many.begin(),expItr=subExport.begin();
00163 vecItr!=many.end();++vecItr,++expItr) {
00164
00165
00166 RCP<const Epetra_MultiVector> srcVec = *vecItr;
00167
00168
00169 const Epetra_BlockMap & globalMap = (*expItr)->SourceMap();
00170
00171
00172 Epetra_MultiVector exportVector(View,globalMap,srcVec->Values(),srcVec->Stride(),srcVec->NumVectors());
00173 one.Export(exportVector,**expItr,Insert);
00174 }
00175 }
00176
00182 RCP<Epetra_IntVector> getSubBlockColumnGIDs(const Epetra_CrsMatrix & A,const MapPair & mapPair)
00183 {
00184 RCP<const Epetra_Map> blkGIDMap = mapPair.first;
00185 RCP<const Epetra_Map> blkContigMap = mapPair.second;
00186
00187
00188 Epetra_IntVector rIndex(A.RowMap(),true);
00189 for(int i=0;i<A.NumMyRows();i++) {
00190
00191 int lid = blkGIDMap->LID(A.GRID(i));
00192 if(lid>-1)
00193 rIndex[i] = blkContigMap->GID(lid);
00194 else
00195 rIndex[i] = -1;
00196 }
00197
00198
00199 Epetra_Import import(A.ColMap(),A.RowMap());
00200 RCP<Epetra_IntVector> cIndex = rcp(new Epetra_IntVector(A.ColMap(),true));
00201 cIndex->Import(rIndex,import,Insert);
00202
00203 return cIndex;
00204 }
00205
00206
00207 RCP<Epetra_CrsMatrix> buildSubBlock(int i,int j,const Epetra_CrsMatrix & A,const std::vector<MapPair> & subMaps)
00208 {
00209
00210 int numVarFamily = subMaps.size();
00211
00212 TEUCHOS_ASSERT(i>=0 && i<numVarFamily);
00213 TEUCHOS_ASSERT(j>=0 && j<numVarFamily);
00214
00215 const Epetra_Map & gRowMap = *subMaps[i].first;
00216 const Epetra_Map & rowMap = *subMaps[i].second;
00217 const Epetra_Map & colMap = *subMaps[j].second;
00218
00219 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A,subMaps[j]);
00220 Epetra_IntVector & local2ContigGIDs = *plocal2ContigGIDs;
00221
00222 RCP<Epetra_CrsMatrix> mat = rcp(new Epetra_CrsMatrix(Copy,rowMap,0));
00223
00224
00225 int numMyRows = rowMap.NumMyElements();
00226 int maxNumEntries = A.MaxNumEntries();
00227
00228
00229 std::vector<int> indices(maxNumEntries);
00230 std::vector<double> values(maxNumEntries);
00231
00232
00233 std::vector<int> colIndices(maxNumEntries);
00234 std::vector<double> colValues(maxNumEntries);
00235
00236
00237
00238 for(int localRow=0;localRow<numMyRows;localRow++) {
00239 int numEntries = -1;
00240 int globalRow = gRowMap.GID(localRow);
00241 int lid = A.RowMap().LID(globalRow);
00242 int contigRow = rowMap.GID(localRow);
00243 TEUCHOS_ASSERT(lid>-1);
00244
00245 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
00246 TEUCHOS_ASSERT(err==0);
00247
00248 int numOwnedCols = 0;
00249 for(int localCol=0;localCol<numEntries;localCol++) {
00250 TEUCHOS_ASSERT(indices[localCol]>-1);
00251
00252
00253 int gid = local2ContigGIDs[indices[localCol]];
00254 if(gid==-1) continue;
00255
00256 colIndices[numOwnedCols] = gid;
00257 colValues[numOwnedCols] = values[localCol];
00258 numOwnedCols++;
00259 }
00260
00261
00262 mat->InsertGlobalValues(contigRow,numOwnedCols,&colValues[0],&colIndices[0]);
00263 }
00264
00265
00266 mat->FillComplete(colMap,rowMap);
00267
00268 return mat;
00269 }
00270
00271
00272 void rebuildSubBlock(int i,int j,const Epetra_CrsMatrix & A,const std::vector<MapPair> & subMaps,Epetra_CrsMatrix & mat)
00273 {
00274
00275 int numVarFamily = subMaps.size();
00276
00277 TEUCHOS_ASSERT(i>=0 && i<numVarFamily);
00278 TEUCHOS_ASSERT(j>=0 && j<numVarFamily);
00279
00280 const Epetra_Map & gRowMap = *subMaps[i].first;
00281 const Epetra_Map & rowMap = *subMaps[i].second;
00282
00283 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A,subMaps[j]);
00284 Epetra_IntVector & local2ContigGIDs = *plocal2ContigGIDs;
00285
00286 mat.PutScalar(0.0);
00287
00288
00289 int numMyRows = rowMap.NumMyElements();
00290 int maxNumEntries = A.MaxNumEntries();
00291
00292
00293 std::vector<int> indices(maxNumEntries);
00294 std::vector<double> values(maxNumEntries);
00295
00296
00297 std::vector<int> colIndices(maxNumEntries);
00298 std::vector<double> colValues(maxNumEntries);
00299
00300
00301
00302 for(int localRow=0;localRow<numMyRows;localRow++) {
00303 int numEntries = -1;
00304 int globalRow = gRowMap.GID(localRow);
00305 int lid = A.RowMap().LID(globalRow);
00306 int contigRow = rowMap.GID(localRow);
00307 TEUCHOS_ASSERT(lid>-1);
00308
00309 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
00310 TEUCHOS_ASSERT(err==0);
00311
00312 int numOwnedCols = 0;
00313 for(int localCol=0;localCol<numEntries;localCol++) {
00314 TEUCHOS_ASSERT(indices[localCol]>-1);
00315
00316
00317 int gid = local2ContigGIDs[indices[localCol]];
00318 if(gid==-1) continue;
00319
00320 colIndices[numOwnedCols] = gid;
00321 colValues[numOwnedCols] = values[localCol];
00322 numOwnedCols++;
00323 }
00324
00325
00326 mat.SumIntoGlobalValues(contigRow,numOwnedCols,&colValues[0],&colIndices[0]);
00327 }
00328 }
00329
00330 }
00331 }
00332 }