47 #include "Teko_BlockingEpetra.hpp"
49 #include "Epetra_IntVector.h"
50 #include "Epetra_LocalMap.h"
72 const MapPair buildSubMap(
const std::vector< int > & gid,
const Epetra_Comm &comm)
74 Teuchos::RCP<Epetra_Map> gidMap = rcp(
new Epetra_Map(-1,gid.size(),&gid[0],0,comm));
75 Teuchos::RCP<Epetra_Map> contigMap = rcp(
new Epetra_Map(-1,gid.size(),0,comm));
77 return std::make_pair(gidMap,contigMap);
89 const ImExPair buildExportImport(
const Epetra_Map & baseMap,
const MapPair & maps)
91 return std::make_pair(rcp(
new Epetra_Import(*maps.first,baseMap)),
92 rcp(
new Epetra_Export(*maps.first,baseMap)));
102 void buildSubVectors(
const std::vector<MapPair> & maps,
103 std::vector<RCP<Epetra_MultiVector> > & vectors,
int count)
105 std::vector<MapPair>::const_iterator mapItr;
108 for(mapItr=maps.begin();mapItr!=maps.end();mapItr++) {
110 RCP<Epetra_MultiVector> mv = rcp(
new Epetra_MultiVector(*(*mapItr).second,count));
111 vectors.push_back(mv);
121 void one2many(std::vector<RCP<Epetra_MultiVector> > & many,
const Epetra_MultiVector & single,
122 const std::vector<RCP<Epetra_Import> > & subImport)
125 std::vector<RCP<Epetra_MultiVector> >::const_iterator vecItr;
126 std::vector<RCP<Epetra_Import> >::const_iterator impItr;
129 for(vecItr=many.begin(),impItr=subImport.begin();
130 vecItr!=many.end();++vecItr,++impItr) {
132 RCP<Epetra_MultiVector> destVec = *vecItr;
135 const Epetra_BlockMap & globalMap = (*impItr)->TargetMap();
138 Epetra_MultiVector importVector(View,globalMap,destVec->Values(),destVec->Stride(),destVec->NumVectors());
141 importVector.Import(single,**impItr,Insert);
154 void many2one(Epetra_MultiVector & one,
const std::vector<RCP<const Epetra_MultiVector> > & many,
155 const std::vector<RCP<Epetra_Export> > & subExport)
158 std::vector<RCP<const Epetra_MultiVector> >::const_iterator vecItr;
159 std::vector<RCP<Epetra_Export> >::const_iterator expItr;
162 for(vecItr=many.begin(),expItr=subExport.begin();
163 vecItr!=many.end();++vecItr,++expItr) {
166 RCP<const Epetra_MultiVector> srcVec = *vecItr;
169 const Epetra_BlockMap & globalMap = (*expItr)->SourceMap();
172 Epetra_MultiVector exportVector(View,globalMap,srcVec->Values(),srcVec->Stride(),srcVec->NumVectors());
173 one.Export(exportVector,**expItr,Insert);
182 RCP<Epetra_IntVector> getSubBlockColumnGIDs(
const Epetra_CrsMatrix & A,
const MapPair & mapPair)
184 RCP<const Epetra_Map> blkGIDMap = mapPair.first;
185 RCP<const Epetra_Map> blkContigMap = mapPair.second;
188 Epetra_IntVector rIndex(A.RowMap(),
true);
189 for(
int i=0;i<A.NumMyRows();i++) {
191 int lid = blkGIDMap->LID(A.GRID(i));
193 rIndex[i] = blkContigMap->GID(lid);
199 Epetra_Import
import(A.ColMap(),A.RowMap());
200 RCP<Epetra_IntVector> cIndex = rcp(
new Epetra_IntVector(A.ColMap(),
true));
201 cIndex->Import(rIndex,
import,Insert);
207 RCP<Epetra_CrsMatrix> buildSubBlock(
int i,
int j,
const Epetra_CrsMatrix & A,
const std::vector<MapPair> & subMaps)
210 int numVarFamily = subMaps.size();
212 TEUCHOS_ASSERT(i>=0 && i<numVarFamily);
213 TEUCHOS_ASSERT(j>=0 && j<numVarFamily);
215 const Epetra_Map & gRowMap = *subMaps[i].first;
216 const Epetra_Map & rowMap = *subMaps[i].second;
217 const Epetra_Map & colMap = *subMaps[j].second;
219 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A,subMaps[j]);
220 Epetra_IntVector & local2ContigGIDs = *plocal2ContigGIDs;
222 RCP<Epetra_CrsMatrix> mat = rcp(
new Epetra_CrsMatrix(Copy,rowMap,0));
225 int numMyRows = rowMap.NumMyElements();
226 int maxNumEntries = A.MaxNumEntries();
229 std::vector<int> indices(maxNumEntries);
230 std::vector<double> values(maxNumEntries);
233 std::vector<int> colIndices(maxNumEntries);
234 std::vector<double> colValues(maxNumEntries);
238 for(
int localRow=0;localRow<numMyRows;localRow++) {
240 int globalRow = gRowMap.GID(localRow);
241 int lid = A.RowMap().LID(globalRow);
242 int contigRow = rowMap.GID(localRow);
243 TEUCHOS_ASSERT(lid>-1);
245 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
246 TEUCHOS_ASSERT(err==0);
248 int numOwnedCols = 0;
249 for(
int localCol=0;localCol<numEntries;localCol++) {
250 TEUCHOS_ASSERT(indices[localCol]>-1);
253 int gid = local2ContigGIDs[indices[localCol]];
254 if(gid==-1)
continue;
256 colIndices[numOwnedCols] = gid;
257 colValues[numOwnedCols] = values[localCol];
262 mat->InsertGlobalValues(contigRow,numOwnedCols,&colValues[0],&colIndices[0]);
266 mat->FillComplete(colMap,rowMap);
272 void rebuildSubBlock(
int i,
int j,
const Epetra_CrsMatrix & A,
const std::vector<MapPair> & subMaps,Epetra_CrsMatrix & mat)
275 int numVarFamily = subMaps.size();
277 TEUCHOS_ASSERT(i>=0 && i<numVarFamily);
278 TEUCHOS_ASSERT(j>=0 && j<numVarFamily);
280 const Epetra_Map & gRowMap = *subMaps[i].first;
281 const Epetra_Map & rowMap = *subMaps[i].second;
283 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A,subMaps[j]);
284 Epetra_IntVector & local2ContigGIDs = *plocal2ContigGIDs;
289 int numMyRows = rowMap.NumMyElements();
290 int maxNumEntries = A.MaxNumEntries();
293 std::vector<int> indices(maxNumEntries);
294 std::vector<double> values(maxNumEntries);
297 std::vector<int> colIndices(maxNumEntries);
298 std::vector<double> colValues(maxNumEntries);
302 for(
int localRow=0;localRow<numMyRows;localRow++) {
304 int globalRow = gRowMap.GID(localRow);
305 int lid = A.RowMap().LID(globalRow);
306 int contigRow = rowMap.GID(localRow);
307 TEUCHOS_ASSERT(lid>-1);
309 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
310 TEUCHOS_ASSERT(err==0);
312 int numOwnedCols = 0;
313 for(
int localCol=0;localCol<numEntries;localCol++) {
314 TEUCHOS_ASSERT(indices[localCol]>-1);
317 int gid = local2ContigGIDs[indices[localCol]];
318 if(gid==-1)
continue;
320 colIndices[numOwnedCols] = gid;
321 colValues[numOwnedCols] = values[localCol];
326 mat.SumIntoGlobalValues(contigRow,numOwnedCols,&colValues[0],&colIndices[0]);