10 #include "Teko_InterlacedEpetra.hpp"
24 void buildSubMaps(
int numGlobals,
int numVars,
const Epetra_Comm& comm,
25 std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps) {
26 std::vector<int> vars;
29 for (
int i = 0; i < numVars; i++) vars.push_back(1);
32 buildSubMaps(numGlobals, vars, comm, subMaps);
36 void buildSubMaps(
const Epetra_Map& globalMap,
const std::vector<int>& vars,
37 const Epetra_Comm& comm,
38 std::vector<std::pair<
int, Teuchos::RCP<Epetra_Map> > >& subMaps) {
39 buildSubMaps(globalMap.NumGlobalElements(), globalMap.NumMyElements(), globalMap.MinMyGID(), vars,
44 void buildSubMaps(
int numGlobals,
const std::vector<int>& vars,
const Epetra_Comm& comm,
45 std::vector<std::pair<
int, Teuchos::RCP<Epetra_Map> > >& subMaps) {
46 std::vector<int>::const_iterator varItr;
49 int numGlobalVars = 0;
50 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) numGlobalVars += *varItr;
53 TEUCHOS_ASSERT((numGlobals % numGlobalVars) == 0);
55 Epetra_Map sampleMap(numGlobals / numGlobalVars, 0, comm);
57 buildSubMaps(numGlobals, numGlobalVars * sampleMap.NumMyElements(),
58 numGlobalVars * sampleMap.MinMyGID(), vars, comm, subMaps);
62 void buildSubMaps(
int numGlobals,
int numMyElements,
int minMyGID,
const std::vector<int>& vars,
63 const Epetra_Comm& comm,
64 std::vector<std::pair<
int, Teuchos::RCP<Epetra_Map> > >& subMaps) {
65 std::vector<int>::const_iterator varItr;
68 int numGlobalVars = 0;
69 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) numGlobalVars += *varItr;
72 TEUCHOS_ASSERT((numGlobals % numGlobalVars) == 0);
73 TEUCHOS_ASSERT((numMyElements % numGlobalVars) == 0);
74 TEUCHOS_ASSERT((minMyGID % numGlobalVars) == 0);
76 int numBlocks = numMyElements / numGlobalVars;
77 int minBlockID = minMyGID / numGlobalVars;
83 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) {
84 int numLocalVars = *varItr;
85 int numAllElmts = numLocalVars * numGlobals / numGlobalVars;
86 int numMyElmts = numLocalVars * numBlocks;
89 std::vector<int> subGlobals;
90 std::vector<int> contigGlobals;
94 for (
int blockNum = 0; blockNum < numBlocks; blockNum++) {
96 for (
int local = 0; local < numLocalVars; ++local) {
100 subGlobals.push_back((minBlockID + blockNum) * numGlobalVars + blockOffset + local);
103 contigGlobals.push_back(numLocalVars * minBlockID + count);
109 assert((
unsigned int)numMyElmts == subGlobals.size());
112 RCP<Epetra_Map> subMap = rcp(
new Epetra_Map(numAllElmts, numMyElmts, &subGlobals[0], 0, comm));
113 RCP<Epetra_Map> contigMap =
114 rcp(
new Epetra_Map(numAllElmts, numMyElmts, &contigGlobals[0], 0, comm));
116 Teuchos::set_extra_data(contigMap,
"contigMap", Teuchos::inOutArg(subMap));
117 subMaps.push_back(std::make_pair(numLocalVars, subMap));
120 blockOffset += numLocalVars;
124 void buildExportImport(
const Epetra_Map& baseMap,
125 const std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps,
126 std::vector<RCP<Epetra_Export> >& subExport,
127 std::vector<RCP<Epetra_Import> >& subImport) {
128 std::vector<std::pair<int, RCP<Epetra_Map> > >::const_iterator mapItr;
131 for (mapItr = subMaps.begin(); mapItr != subMaps.end(); ++mapItr) {
133 const Epetra_Map& map = *(mapItr->second);
136 subImport.push_back(rcp(
new Epetra_Import(map, baseMap)));
137 subExport.push_back(rcp(
new Epetra_Export(map, baseMap)));
141 void buildSubVectors(
const std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps,
142 std::vector<RCP<Epetra_MultiVector> >& subVectors,
int count) {
143 std::vector<std::pair<int, RCP<Epetra_Map> > >::const_iterator mapItr;
146 for (mapItr = subMaps.begin(); mapItr != subMaps.end(); ++mapItr) {
148 const Epetra_Map& map =
149 *(Teuchos::get_extra_data<RCP<Epetra_Map> >(mapItr->second,
"contigMap"));
152 RCP<Epetra_MultiVector> mv = rcp(
new Epetra_MultiVector(map, count));
153 Teuchos::set_extra_data(mapItr->second,
"globalMap", Teuchos::inOutArg(mv));
154 subVectors.push_back(mv);
158 void associateSubVectors(
const std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps,
159 std::vector<RCP<const Epetra_MultiVector> >& subVectors) {
160 std::vector<std::pair<int, RCP<Epetra_Map> > >::const_iterator mapItr;
161 std::vector<RCP<const Epetra_MultiVector> >::iterator vecItr;
163 TEUCHOS_ASSERT(subMaps.size() == subVectors.size());
166 for (mapItr = subMaps.begin(), vecItr = subVectors.begin(); mapItr != subMaps.end();
168 Teuchos::set_extra_data(mapItr->second,
"globalMap", Teuchos::inOutArg(*vecItr));
172 RCP<Epetra_CrsMatrix> buildSubBlock(
int i,
int j,
const Epetra_CrsMatrix& A,
173 const std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps) {
175 int numVarFamily = subMaps.size();
177 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
178 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
180 const Epetra_Map& gRowMap = *subMaps[i].second;
181 const Epetra_Map& rowMap =
182 *Teuchos::get_extra_data<RCP<Epetra_Map> >(subMaps[i].second,
"contigMap");
183 const Epetra_Map& colMap =
184 *Teuchos::get_extra_data<RCP<Epetra_Map> >(subMaps[j].second,
"contigMap");
185 int colFamilyCnt = subMaps[j].first;
189 int numGlobalVars = 0;
190 int rowBlockOffset = 0;
191 int colBlockOffset = 0;
192 for (
int k = 0; k < numVarFamily; k++) {
193 numGlobalVars += subMaps[k].first;
196 if (k < i) rowBlockOffset += subMaps[k].first;
197 if (k < j) colBlockOffset += subMaps[k].first;
201 Epetra_Import
import(gRowMap, A.RowMap());
202 Epetra_CrsMatrix localA(Copy, gRowMap, 0);
203 localA.Import(A,
import, Insert);
205 RCP<Epetra_CrsMatrix> mat = rcp(
new Epetra_CrsMatrix(Copy, rowMap, 0));
208 int numMyRows = rowMap.NumMyElements();
209 int maxNumEntries = A.GlobalMaxNumEntries();
212 std::vector<int> indices(maxNumEntries);
213 std::vector<double> values(maxNumEntries);
216 std::vector<int> colIndices(maxNumEntries);
217 std::vector<double> colValues(maxNumEntries);
221 for (
int localRow = 0; localRow < numMyRows; localRow++) {
223 int globalRow = gRowMap.GID(localRow);
224 int contigRow = rowMap.GID(localRow);
226 TEUCHOS_ASSERT(globalRow >= 0);
227 TEUCHOS_ASSERT(contigRow >= 0);
231 localA.ExtractGlobalRowCopy(globalRow, maxNumEntries, numEntries, &values[0], &indices[0]);
232 TEUCHOS_ASSERT(err == 0);
234 int numOwnedCols = 0;
235 for (
int localCol = 0; localCol < numEntries; localCol++) {
236 int globalCol = indices[localCol];
239 int block = globalCol / numGlobalVars;
241 bool inFamily =
true;
244 inFamily &= (block * numGlobalVars + colBlockOffset <= globalCol);
245 inFamily &= ((block * numGlobalVars + colBlockOffset + colFamilyCnt) > globalCol);
249 int familyOffset = globalCol - (block * numGlobalVars + colBlockOffset);
251 colIndices[numOwnedCols] = block * colFamilyCnt + familyOffset;
252 colValues[numOwnedCols] = values[localCol];
259 mat->InsertGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
263 mat->FillComplete(colMap, rowMap);
269 void rebuildSubBlock(
int i,
int j,
const Epetra_CrsMatrix& A,
270 const std::vector<std::pair<
int, RCP<Epetra_Map> > >& subMaps,
271 Epetra_CrsMatrix& mat) {
273 int numVarFamily = subMaps.size();
275 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
276 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
277 TEUCHOS_ASSERT(mat.Filled());
279 const Epetra_Map& gRowMap = *subMaps[i].second;
280 const Epetra_Map& rowMap =
281 *Teuchos::get_extra_data<RCP<Epetra_Map> >(subMaps[i].second,
"contigMap");
282 int colFamilyCnt = subMaps[j].first;
286 int numGlobalVars = 0;
287 int rowBlockOffset = 0;
288 int colBlockOffset = 0;
289 for (
int k = 0; k < numVarFamily; k++) {
290 numGlobalVars += subMaps[k].first;
293 if (k < i) rowBlockOffset += subMaps[k].first;
294 if (k < j) colBlockOffset += subMaps[k].first;
298 Epetra_Import
import(gRowMap, A.RowMap());
299 Epetra_CrsMatrix localA(Copy, gRowMap, 0);
300 localA.Import(A,
import, Insert);
306 int numMyRows = rowMap.NumMyElements();
307 int maxNumEntries = A.GlobalMaxNumEntries();
310 std::vector<int> indices(maxNumEntries);
311 std::vector<double> values(maxNumEntries);
314 std::vector<int> colIndices(maxNumEntries);
315 std::vector<double> colValues(maxNumEntries);
319 for (
int localRow = 0; localRow < numMyRows; localRow++) {
321 int globalRow = gRowMap.GID(localRow);
322 int contigRow = rowMap.GID(localRow);
324 TEUCHOS_ASSERT(globalRow >= 0);
325 TEUCHOS_ASSERT(contigRow >= 0);
329 localA.ExtractGlobalRowCopy(globalRow, maxNumEntries, numEntries, &values[0], &indices[0]);
330 TEUCHOS_ASSERT(err == 0);
332 int numOwnedCols = 0;
333 for (
int localCol = 0; localCol < numEntries; localCol++) {
334 int globalCol = indices[localCol];
337 int block = globalCol / numGlobalVars;
339 bool inFamily =
true;
342 inFamily &= (block * numGlobalVars + colBlockOffset <= globalCol);
343 inFamily &= ((block * numGlobalVars + colBlockOffset + colFamilyCnt) > globalCol);
347 int familyOffset = globalCol - (block * numGlobalVars + colBlockOffset);
349 colIndices[numOwnedCols] = block * colFamilyCnt + familyOffset;
350 colValues[numOwnedCols] = values[localCol];
357 mat.SumIntoGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
362 void many2one(Epetra_MultiVector& one,
const std::vector<RCP<const Epetra_MultiVector> >& many,
363 const std::vector<RCP<Epetra_Export> >& subExport) {
365 std::vector<RCP<const Epetra_MultiVector> >::const_iterator vecItr;
366 std::vector<RCP<Epetra_Export> >::const_iterator expItr;
369 for (vecItr = many.begin(), expItr = subExport.begin(); vecItr != many.end();
370 ++vecItr, ++expItr) {
372 RCP<const Epetra_MultiVector> srcVec = *vecItr;
375 const Epetra_Map& globalMap = *(Teuchos::get_extra_data<RCP<Epetra_Map> >(srcVec,
"globalMap"));
378 Epetra_MultiVector exportVector(View, globalMap, srcVec->Values(), srcVec->Stride(),
379 srcVec->NumVectors());
380 one.Export(exportVector, **expItr, Insert);
385 void one2many(std::vector<RCP<Epetra_MultiVector> >& many,
const Epetra_MultiVector& single,
386 const std::vector<RCP<Epetra_Import> >& subImport) {
388 std::vector<RCP<Epetra_MultiVector> >::const_iterator vecItr;
389 std::vector<RCP<Epetra_Import> >::const_iterator impItr;
392 for (vecItr = many.begin(), impItr = subImport.begin(); vecItr != many.end();
393 ++vecItr, ++impItr) {
395 RCP<Epetra_MultiVector> destVec = *vecItr;
398 const Epetra_Map& globalMap =
399 *(Teuchos::get_extra_data<RCP<Epetra_Map> >(destVec,
"globalMap"));
402 Epetra_MultiVector importVector(View, globalMap, destVec->Values(), destVec->Stride(),
403 destVec->NumVectors());
406 importVector.Import(single, **impItr, Insert);