Teko Version of the Day
Loading...
Searching...
No Matches
Teko_BlockingEpetra.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_BlockingEpetra.hpp"
11
12#include "Epetra_IntVector.h"
13#include "Epetra_LocalMap.h"
14
15using Teuchos::RCP;
16using Teuchos::rcp;
17
18namespace Teko {
19namespace Epetra {
20namespace Blocking {
21
35const MapPair buildSubMap(const std::vector<int> &gid, const Epetra_Comm &comm) {
36 Teuchos::RCP<Epetra_Map> gidMap = rcp(new Epetra_Map(-1, gid.size(), &gid[0], 0, comm));
37 Teuchos::RCP<Epetra_Map> contigMap = rcp(new Epetra_Map(-1, gid.size(), 0, comm));
38
39 return std::make_pair(gidMap, contigMap);
40}
41
51const ImExPair buildExportImport(const Epetra_Map &baseMap, const MapPair &maps) {
52 return std::make_pair(rcp(new Epetra_Import(*maps.first, baseMap)),
53 rcp(new Epetra_Export(*maps.first, baseMap)));
54}
55
63void buildSubVectors(const std::vector<MapPair> &maps,
64 std::vector<RCP<Epetra_MultiVector> > &vectors, int count) {
65 std::vector<MapPair>::const_iterator mapItr;
66
67 // loop over all maps
68 for (mapItr = maps.begin(); mapItr != maps.end(); mapItr++) {
69 // add new elements to vectors
70 RCP<Epetra_MultiVector> mv = rcp(new Epetra_MultiVector(*(*mapItr).second, count));
71 vectors.push_back(mv);
72 }
73}
74
82void one2many(std::vector<RCP<Epetra_MultiVector> > &many, const Epetra_MultiVector &single,
83 const std::vector<RCP<Epetra_Import> > &subImport) {
84 // std::vector<RCP<Epetra_Vector> >::const_iterator vecItr;
85 std::vector<RCP<Epetra_MultiVector> >::const_iterator vecItr;
86 std::vector<RCP<Epetra_Import> >::const_iterator impItr;
87
88 // using Importers fill the sub vectors from the mama vector
89 for (vecItr = many.begin(), impItr = subImport.begin(); vecItr != many.end();
90 ++vecItr, ++impItr) {
91 // for ease of access to the destination
92 RCP<Epetra_MultiVector> destVec = *vecItr;
93
94 // extract the map with global indicies from the current vector
95 const Epetra_BlockMap &globalMap = (*impItr)->TargetMap();
96
97 // build the import vector as a view on the destination
98 Epetra_MultiVector importVector(View, globalMap, destVec->Values(), destVec->Stride(),
99 destVec->NumVectors());
100
101 // perform the import
102 importVector.Import(single, **impItr, Insert);
103 }
104}
105
116void many2one(Epetra_MultiVector &one, const std::vector<RCP<const Epetra_MultiVector> > &many,
117 const std::vector<RCP<Epetra_Export> > &subExport) {
118 // std::vector<RCP<const Epetra_Vector> >::const_iterator vecItr;
119 std::vector<RCP<const Epetra_MultiVector> >::const_iterator vecItr;
120 std::vector<RCP<Epetra_Export> >::const_iterator expItr;
121
122 // using Exporters fill the empty vector from the sub-vectors
123 for (vecItr = many.begin(), expItr = subExport.begin(); vecItr != many.end();
124 ++vecItr, ++expItr) {
125 // for ease of access to the source
126 RCP<const Epetra_MultiVector> srcVec = *vecItr;
127
128 // extract the map with global indicies from the current vector
129 const Epetra_BlockMap &globalMap = (*expItr)->SourceMap();
130
131 // build the export vector as a view of the destination
132 Epetra_MultiVector exportVector(View, globalMap, srcVec->Values(), srcVec->Stride(),
133 srcVec->NumVectors());
134 one.Export(exportVector, **expItr, Insert);
135 }
136}
137
143RCP<Epetra_IntVector> getSubBlockColumnGIDs(const Epetra_CrsMatrix &A, const MapPair &mapPair) {
144 RCP<const Epetra_Map> blkGIDMap = mapPair.first;
145 RCP<const Epetra_Map> blkContigMap = mapPair.second;
146
147 // fill index vector for rows
148 Epetra_IntVector rIndex(A.RowMap(), true);
149 for (int i = 0; i < A.NumMyRows(); i++) {
150 // LID is need to get to contiguous map
151 int lid = blkGIDMap->LID(A.GRID(i)); // this LID makes me nervous
152 if (lid > -1)
153 rIndex[i] = blkContigMap->GID(lid);
154 else
155 rIndex[i] = -1;
156 }
157
158 // get relavant column indices
159 Epetra_Import import(A.ColMap(), A.RowMap());
160 RCP<Epetra_IntVector> cIndex = rcp(new Epetra_IntVector(A.ColMap(), true));
161 cIndex->Import(rIndex, import, Insert);
162
163 return cIndex;
164}
165
166// build a single subblock Epetra_CrsMatrix
167RCP<Epetra_CrsMatrix> buildSubBlock(int i, int j, const Epetra_CrsMatrix &A,
168 const std::vector<MapPair> &subMaps) {
169 // get the number of variables families
170 int numVarFamily = subMaps.size();
171
172 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
173 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
174
175 const Epetra_Map &gRowMap = *subMaps[i].first; // new GIDs
176 const Epetra_Map &rowMap = *subMaps[i].second; // contiguous GIDs
177 const Epetra_Map &colMap = *subMaps[j].second;
178
179 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A, subMaps[j]);
180 Epetra_IntVector &local2ContigGIDs = *plocal2ContigGIDs;
181
182 RCP<Epetra_CrsMatrix> mat = rcp(new Epetra_CrsMatrix(Copy, rowMap, 0));
183
184 // get entry information
185 int numMyRows = rowMap.NumMyElements();
186 int maxNumEntries = A.MaxNumEntries();
187
188 // for extraction
189 std::vector<int> indices(maxNumEntries);
190 std::vector<double> values(maxNumEntries);
191
192 // for insertion
193 std::vector<int> colIndices(maxNumEntries);
194 std::vector<double> colValues(maxNumEntries);
195
196 // insert each row into subblock
197 // let FillComplete handle column distribution
198 for (int localRow = 0; localRow < numMyRows; localRow++) {
199 int numEntries = -1;
200 int globalRow = gRowMap.GID(localRow);
201 int lid = A.RowMap().LID(globalRow);
202 int contigRow = rowMap.GID(localRow);
203 TEUCHOS_ASSERT(lid > -1);
204
205 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
206 TEUCHOS_ASSERT(err == 0);
207
208 int numOwnedCols = 0;
209 for (int localCol = 0; localCol < numEntries; localCol++) {
210 TEUCHOS_ASSERT(indices[localCol] > -1);
211
212 // if global id is not owned by this column
213 int gid = local2ContigGIDs[indices[localCol]];
214 if (gid == -1) continue; // in contiguous row
215
216 colIndices[numOwnedCols] = gid;
217 colValues[numOwnedCols] = values[localCol];
218 numOwnedCols++;
219 }
220
221 // insert it into the new matrix
222 mat->InsertGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
223 }
224
225 // fill it and automagically optimize the storage
226 mat->FillComplete(colMap, rowMap);
227
228 return mat;
229}
230
231// build a single subblock Epetra_CrsMatrix
232void rebuildSubBlock(int i, int j, const Epetra_CrsMatrix &A, const std::vector<MapPair> &subMaps,
233 Epetra_CrsMatrix &mat) {
234 // get the number of variables families
235 int numVarFamily = subMaps.size();
236
237 TEUCHOS_ASSERT(i >= 0 && i < numVarFamily);
238 TEUCHOS_ASSERT(j >= 0 && j < numVarFamily);
239
240 const Epetra_Map &gRowMap = *subMaps[i].first; // new GIDs
241 const Epetra_Map &rowMap = *subMaps[i].second; // contiguous GIDs
242
243 const RCP<Epetra_IntVector> plocal2ContigGIDs = getSubBlockColumnGIDs(A, subMaps[j]);
244 Epetra_IntVector &local2ContigGIDs = *plocal2ContigGIDs;
245
246 mat.PutScalar(0.0);
247
248 // get entry information
249 int numMyRows = rowMap.NumMyElements();
250 int maxNumEntries = A.MaxNumEntries();
251
252 // for extraction
253 std::vector<int> indices(maxNumEntries);
254 std::vector<double> values(maxNumEntries);
255
256 // for insertion
257 std::vector<int> colIndices(maxNumEntries);
258 std::vector<double> colValues(maxNumEntries);
259
260 // insert each row into subblock
261 // let FillComplete handle column distribution
262 for (int localRow = 0; localRow < numMyRows; localRow++) {
263 int numEntries = -1;
264 int globalRow = gRowMap.GID(localRow);
265 int lid = A.RowMap().LID(globalRow);
266 int contigRow = rowMap.GID(localRow);
267 TEUCHOS_ASSERT(lid > -1);
268
269 int err = A.ExtractMyRowCopy(lid, maxNumEntries, numEntries, &values[0], &indices[0]);
270 TEUCHOS_ASSERT(err == 0);
271
272 int numOwnedCols = 0;
273 for (int localCol = 0; localCol < numEntries; localCol++) {
274 TEUCHOS_ASSERT(indices[localCol] > -1);
275
276 // if global id is not owned by this column
277 int gid = local2ContigGIDs[indices[localCol]];
278 if (gid == -1) continue; // in contiguous row
279
280 colIndices[numOwnedCols] = gid;
281 colValues[numOwnedCols] = values[localCol];
282 numOwnedCols++;
283 }
284
285 // insert it into the new matrix
286 mat.SumIntoGlobalValues(contigRow, numOwnedCols, &colValues[0], &colIndices[0]);
287 }
288}
289
290} // namespace Blocking
291} // namespace Epetra
292} // namespace Teko