Teko Version of the Day
Loading...
Searching...
No Matches
Teko_InterlacedEpetra.cpp
1// @HEADER
2// *****************************************************************************
3// Teko: A package for block and physics based preconditioning
4//
5// Copyright 2010 NTESS and the Teko contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#include "Teko_InterlacedEpetra.hpp"
11
12#include <vector>
13
14using Teuchos::RCP;
15using Teuchos::rcp;
16
17namespace Teko {
18namespace Epetra {
19namespace Strided {
20
21// this assumes that there are numGlobals with numVars each interlaced
22// i.e. for numVars = 2 (u,v) then the vector is
23// [u_0,v_0,u_1,v_1,u_2,v_2, ..., u_(numGlobals-1),v_(numGlobals-1)]
24void buildSubMaps(int numGlobals, int numVars, const Epetra_Comm& comm,
25 std::vector<std::pair<int, RCP<Epetra_Map> > >& subMaps) {
26 std::vector<int> vars;
27
28 // build vector describing the sub maps
29 for (int i = 0; i < numVars; i++) vars.push_back(1);
30
31 // build all the submaps
32 buildSubMaps(numGlobals, vars, comm, subMaps);
33}
34
35// build maps to make other conversions
36void 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,
40 comm, subMaps);
41}
42
43// build maps to make other conversions
44void 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;
47
48 // compute total number of variables
49 int numGlobalVars = 0;
50 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) numGlobalVars += *varItr;
51
52 // must be an even number of globals
53 TEUCHOS_ASSERT((numGlobals % numGlobalVars) == 0);
54
55 Epetra_Map sampleMap(numGlobals / numGlobalVars, 0, comm);
56
57 buildSubMaps(numGlobals, numGlobalVars * sampleMap.NumMyElements(),
58 numGlobalVars * sampleMap.MinMyGID(), vars, comm, subMaps);
59}
60
61// build maps to make other conversions
62void 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;
66
67 // compute total number of variables
68 int numGlobalVars = 0;
69 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) numGlobalVars += *varItr;
70
71 // must be an even number of globals
72 TEUCHOS_ASSERT((numGlobals % numGlobalVars) == 0);
73 TEUCHOS_ASSERT((numMyElements % numGlobalVars) == 0);
74 TEUCHOS_ASSERT((minMyGID % numGlobalVars) == 0);
75
76 int numBlocks = numMyElements / numGlobalVars;
77 int minBlockID = minMyGID / numGlobalVars;
78
79 subMaps.clear();
80
81 // index into local block in strided map
82 int blockOffset = 0;
83 for (varItr = vars.begin(); varItr != vars.end(); ++varItr) {
84 int numLocalVars = *varItr;
85 int numAllElmts = numLocalVars * numGlobals / numGlobalVars;
86 int numMyElmts = numLocalVars * numBlocks;
87
88 // create global arrays describing the as of yet uncreated maps
89 std::vector<int> subGlobals;
90 std::vector<int> contigGlobals; // the contiguous globals
91
92 // loop over each block of variables
93 int count = 0;
94 for (int blockNum = 0; blockNum < numBlocks; blockNum++) {
95 // loop over each local variable in the block
96 for (int local = 0; local < numLocalVars; ++local) {
97 // global block number = minGID+blockNum
98 // block begin global id = numGlobalVars*(minGID+blockNum)
99 // global id block offset = blockOffset+local
100 subGlobals.push_back((minBlockID + blockNum) * numGlobalVars + blockOffset + local);
101
102 // also build the contiguous IDs
103 contigGlobals.push_back(numLocalVars * minBlockID + count);
104 count++;
105 }
106 }
107
108 // sanity check
109 assert((unsigned int)numMyElmts == subGlobals.size());
110
111 // create the map with contiguous elements and the map with global elements
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));
115
116 Teuchos::set_extra_data(contigMap, "contigMap", Teuchos::inOutArg(subMap));
117 subMaps.push_back(std::make_pair(numLocalVars, subMap));
118
119 // update the block offset
120 blockOffset += numLocalVars;
121 }
122}
123
124void 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;
129
130 // build importers and exporters
131 for (mapItr = subMaps.begin(); mapItr != subMaps.end(); ++mapItr) {
132 // exctract basic map
133 const Epetra_Map& map = *(mapItr->second);
134
135 // add new elements to vectors
136 subImport.push_back(rcp(new Epetra_Import(map, baseMap)));
137 subExport.push_back(rcp(new Epetra_Export(map, baseMap)));
138 }
139}
140
141void 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;
144
145 // build vectors
146 for (mapItr = subMaps.begin(); mapItr != subMaps.end(); ++mapItr) {
147 // exctract basic map
148 const Epetra_Map& map =
149 *(Teuchos::get_extra_data<RCP<Epetra_Map> >(mapItr->second, "contigMap"));
150
151 // add new elements to vectors
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);
155 }
156}
157
158void 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;
162
163 TEUCHOS_ASSERT(subMaps.size() == subVectors.size());
164
165 // associate the sub vectors with the subMaps
166 for (mapItr = subMaps.begin(), vecItr = subVectors.begin(); mapItr != subMaps.end();
167 ++mapItr, ++vecItr)
168 Teuchos::set_extra_data(mapItr->second, "globalMap", Teuchos::inOutArg(*vecItr));
169}
170
171// build a single subblock Epetra_CrsMatrix
172RCP<Epetra_CrsMatrix> buildSubBlock(int i, int j, const Epetra_CrsMatrix& A,
173 const std::vector<std::pair<int, RCP<Epetra_Map> > >& subMaps) {
174 // get the number of variables families
175 int numVarFamily = subMaps.size();
176
177 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
178 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
179
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;
186
187 // compute the number of global variables
188 // and the row and column block offset
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;
194
195 // compute block offsets
196 if (k < i) rowBlockOffset += subMaps[k].first;
197 if (k < j) colBlockOffset += subMaps[k].first;
198 }
199
200 // copy all global rows to here
201 Epetra_Import import(gRowMap, A.RowMap());
202 Epetra_CrsMatrix localA(Copy, gRowMap, 0);
203 localA.Import(A, import, Insert);
204
205 RCP<Epetra_CrsMatrix> mat = rcp(new Epetra_CrsMatrix(Copy, rowMap, 0));
206
207 // get entry information
208 int numMyRows = rowMap.NumMyElements();
209 int maxNumEntries = A.GlobalMaxNumEntries();
210
211 // for extraction
212 std::vector<int> indices(maxNumEntries);
213 std::vector<double> values(maxNumEntries);
214
215 // for insertion
216 std::vector<int> colIndices(maxNumEntries);
217 std::vector<double> colValues(maxNumEntries);
218
219 // insert each row into subblock
220 // let FillComplete handle column distribution
221 for (int localRow = 0; localRow < numMyRows; localRow++) {
222 int numEntries = -1;
223 int globalRow = gRowMap.GID(localRow);
224 int contigRow = rowMap.GID(localRow);
225
226 TEUCHOS_ASSERT(globalRow >= 0);
227 TEUCHOS_ASSERT(contigRow >= 0);
228
229 // extract a global row copy
230 int err =
231 localA.ExtractGlobalRowCopy(globalRow, maxNumEntries, numEntries, &values[0], &indices[0]);
232 TEUCHOS_ASSERT(err == 0);
233
234 int numOwnedCols = 0;
235 for (int localCol = 0; localCol < numEntries; localCol++) {
236 int globalCol = indices[localCol];
237
238 // determinate which block this column ID is in
239 int block = globalCol / numGlobalVars;
240
241 bool inFamily = true;
242
243 // test the beginning of the block
244 inFamily &= (block * numGlobalVars + colBlockOffset <= globalCol);
245 inFamily &= ((block * numGlobalVars + colBlockOffset + colFamilyCnt) > globalCol);
246
247 // is this column in the variable family
248 if (inFamily) {
249 int familyOffset = globalCol - (block * numGlobalVars + colBlockOffset);
250
251 colIndices[numOwnedCols] = block * colFamilyCnt + familyOffset;
252 colValues[numOwnedCols] = values[localCol];
253
254 numOwnedCols++;
255 }
256 }
257
258 // insert it into the new matrix
259 mat->InsertGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
260 }
261
262 // fill it and automagically optimize the storage
263 mat->FillComplete(colMap, rowMap);
264
265 return mat;
266}
267
268// rebuild a single subblock Epetra_CrsMatrix
269void rebuildSubBlock(int i, int j, const Epetra_CrsMatrix& A,
270 const std::vector<std::pair<int, RCP<Epetra_Map> > >& subMaps,
271 Epetra_CrsMatrix& mat) {
272 // get the number of variables families
273 int numVarFamily = subMaps.size();
274
275 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
276 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
277 TEUCHOS_ASSERT(mat.Filled());
278
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;
283
284 // compute the number of global variables
285 // and the row and column block offset
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;
291
292 // compute block offsets
293 if (k < i) rowBlockOffset += subMaps[k].first;
294 if (k < j) colBlockOffset += subMaps[k].first;
295 }
296
297 // copy all global rows to here
298 Epetra_Import import(gRowMap, A.RowMap());
299 Epetra_CrsMatrix localA(Copy, gRowMap, 0);
300 localA.Import(A, import, Insert);
301
302 // clear out the old matrix
303 mat.PutScalar(0.0);
304
305 // get entry information
306 int numMyRows = rowMap.NumMyElements();
307 int maxNumEntries = A.GlobalMaxNumEntries();
308
309 // for extraction
310 std::vector<int> indices(maxNumEntries);
311 std::vector<double> values(maxNumEntries);
312
313 // for insertion
314 std::vector<int> colIndices(maxNumEntries);
315 std::vector<double> colValues(maxNumEntries);
316
317 // insert each row into subblock
318 // let FillComplete handle column distribution
319 for (int localRow = 0; localRow < numMyRows; localRow++) {
320 int numEntries = -1;
321 int globalRow = gRowMap.GID(localRow);
322 int contigRow = rowMap.GID(localRow);
323
324 TEUCHOS_ASSERT(globalRow >= 0);
325 TEUCHOS_ASSERT(contigRow >= 0);
326
327 // extract a global row copy
328 int err =
329 localA.ExtractGlobalRowCopy(globalRow, maxNumEntries, numEntries, &values[0], &indices[0]);
330 TEUCHOS_ASSERT(err == 0);
331
332 int numOwnedCols = 0;
333 for (int localCol = 0; localCol < numEntries; localCol++) {
334 int globalCol = indices[localCol];
335
336 // determinate which block this column ID is in
337 int block = globalCol / numGlobalVars;
338
339 bool inFamily = true;
340
341 // test the beginning of the block
342 inFamily &= (block * numGlobalVars + colBlockOffset <= globalCol);
343 inFamily &= ((block * numGlobalVars + colBlockOffset + colFamilyCnt) > globalCol);
344
345 // is this column in the variable family
346 if (inFamily) {
347 int familyOffset = globalCol - (block * numGlobalVars + colBlockOffset);
348
349 colIndices[numOwnedCols] = block * colFamilyCnt + familyOffset;
350 colValues[numOwnedCols] = values[localCol];
351
352 numOwnedCols++;
353 }
354 }
355
356 // insert it into the new matrix
357 mat.SumIntoGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
358 }
359}
360
361// collect subvectors into a single global vector
362void many2one(Epetra_MultiVector& one, const std::vector<RCP<const Epetra_MultiVector> >& many,
363 const std::vector<RCP<Epetra_Export> >& subExport) {
364 // std::vector<RCP<const Epetra_Vector> >::const_iterator vecItr;
365 std::vector<RCP<const Epetra_MultiVector> >::const_iterator vecItr;
366 std::vector<RCP<Epetra_Export> >::const_iterator expItr;
367
368 // using Exporters fill the empty vector from the sub-vectors
369 for (vecItr = many.begin(), expItr = subExport.begin(); vecItr != many.end();
370 ++vecItr, ++expItr) {
371 // for ease of access to the source
372 RCP<const Epetra_MultiVector> srcVec = *vecItr;
373
374 // extract the map with global indicies from the current vector
375 const Epetra_Map& globalMap = *(Teuchos::get_extra_data<RCP<Epetra_Map> >(srcVec, "globalMap"));
376
377 // build the export vector as a view of the destination
378 Epetra_MultiVector exportVector(View, globalMap, srcVec->Values(), srcVec->Stride(),
379 srcVec->NumVectors());
380 one.Export(exportVector, **expItr, Insert);
381 }
382}
383
384// distribute one global vector into a many subvectors
385void one2many(std::vector<RCP<Epetra_MultiVector> >& many, const Epetra_MultiVector& single,
386 const std::vector<RCP<Epetra_Import> >& subImport) {
387 // std::vector<RCP<Epetra_Vector> >::const_iterator vecItr;
388 std::vector<RCP<Epetra_MultiVector> >::const_iterator vecItr;
389 std::vector<RCP<Epetra_Import> >::const_iterator impItr;
390
391 // using Importers fill the sub vectors from the mama vector
392 for (vecItr = many.begin(), impItr = subImport.begin(); vecItr != many.end();
393 ++vecItr, ++impItr) {
394 // for ease of access to the destination
395 RCP<Epetra_MultiVector> destVec = *vecItr;
396
397 // extract the map with global indicies from the current vector
398 const Epetra_Map& globalMap =
399 *(Teuchos::get_extra_data<RCP<Epetra_Map> >(destVec, "globalMap"));
400
401 // build the import vector as a view on the destination
402 Epetra_MultiVector importVector(View, globalMap, destVec->Values(), destVec->Stride(),
403 destVec->NumVectors());
404
405 // perform the import
406 importVector.Import(single, **impItr, Insert);
407 }
408}
409
410} // namespace Strided
411} // end namespace Epetra
412} // end namespace Teko