Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_getNumDiags.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_DETAILS_GETNUMDIAGS_HPP
11#define TPETRA_DETAILS_GETNUMDIAGS_HPP
12
19
20#include "Tpetra_CrsGraph.hpp"
21#include "Teuchos_CommHelpers.hpp"
23
24namespace Tpetra {
25namespace Details {
26
27namespace Impl {
33 template<class LocalGraphType, class LocalMapType>
34 class CountLocalNumDiags {
35 public:
36 CountLocalNumDiags (const LocalGraphType& G,
37 const LocalMapType& rowMap,
38 const LocalMapType& colMap) :
39 G_ (G), rowMap_ (rowMap), colMap_ (colMap)
40 {}
41
42 // Result can't be more than the number of local rows, so
43 // local_ordinal_type is appropriate.
44 using result_type = typename LocalMapType::local_ordinal_type;
45
47 KOKKOS_INLINE_FUNCTION void
48 operator () (const typename LocalMapType::local_ordinal_type lclRow,
49 result_type& diagCount) const
50 {
51 using LO = typename LocalMapType::local_ordinal_type;
52 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
53
54 auto G_row = G_.rowConst (lclRow);
55 const LO numEnt = G_row.length;
56 if (numEnt != 0) {
57 // Use global row and column indices to find the diagonal
58 // entry. Caller promises that local row index is in the row
59 // Map on the calling process.
60 const LO lclDiagCol = colMap_.getLocalElement (rowMap_.getGlobalElement (lclRow));
61 // If it's not in the column Map, then there's no diagonal entry.
62 if (lclDiagCol != LOT::invalid ()) {
63 // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize for
64 // the sorted case, but note that it requires operator[].
65 for (LO k = 0; k < numEnt; ++k) {
66 if (lclDiagCol == G_row(k)) {
67 ++diagCount;
68 break; // don't count duplicates
69 }
70 }
71 }
72 }
73 }
74
75 private:
76 LocalGraphType G_;
77 LocalMapType rowMap_;
78 LocalMapType colMap_;
79 };
80
81 template<class LO, class GO, class NT>
82 typename ::Tpetra::CrsGraph<LO, GO, NT>::local_ordinal_type
83 countLocalNumDiagsInFillCompleteGraph (const ::Tpetra::CrsGraph<LO, GO, NT>& G)
84 {
86 using local_map_type = typename crs_graph_type::map_type::local_map_type;
87 using local_graph_device_type = typename crs_graph_type::local_graph_device_type;
88 using functor_type = CountLocalNumDiags<local_graph_device_type, local_map_type>;
89 using execution_space = typename crs_graph_type::device_type::execution_space;
90 using policy_type = Kokkos::RangePolicy<execution_space, LO>;
91
92 const auto rowMap = G.getRowMap ();
93 const auto colMap = G.getColMap ();
94 if (rowMap.get () == nullptr || colMap.get () == nullptr) {
95 return 0; // this process does not participate
96 }
97 else {
98 LO lclNumDiags {0};
99 functor_type f (G.getLocalGraphDevice (), rowMap->getLocalMap (), colMap->getLocalMap ());
100 Kokkos::parallel_reduce (policy_type (0, G.getLocalNumRows ()), f, lclNumDiags);
101 return lclNumDiags;
102 }
103 }
104
114 template<class MapType>
115 typename MapType::local_ordinal_type
116 getLocalDiagonalColumnIndex (const typename MapType::local_ordinal_type lclRow,
117 const MapType& rowMap,
118 const MapType& colMap)
119 {
120 return colMap.getLocalElement (rowMap.getGlobalElement (lclRow));
121 }
122
124 template<class LO, class GO, class NT>
125 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
127 {
128 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
129
130 const auto rowMap = G.getRowMap ();
131 const auto colMap = G.getColMap ();
132 if (rowMap.get () == nullptr || colMap.get () == nullptr) {
133 return 0; // this process does not participate
134 }
135 else {
136 TEUCHOS_TEST_FOR_EXCEPTION
137 (! G.supportsRowViews (), std::logic_error, "Not implemented!");
138
139 typename ::Tpetra::RowGraph<LO, GO, NT>::local_inds_host_view_type
140 lclColInds;
141 const LO lclNumRows = static_cast<LO> (G.getLocalNumRows ());
142
143 LO diagCount = 0;
144 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
145 G.getLocalRowView (lclRow, lclColInds);
146 const LO numEnt = static_cast<LO> (lclColInds.size ());
147 if (numEnt != 0) {
148 const LO lclDiagCol = colMap->getLocalElement (rowMap->getGlobalElement (lclRow));
149 // If it's not in the column Map, then there's no diagonal entry.
150 if (lclDiagCol != LOT::invalid ()) {
151 // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize
152 // for the sorted case.
153 for (LO k = 0; k < numEnt; ++k) {
154 if (lclDiagCol == lclColInds[k]) {
155 ++diagCount;
156 break; // don't count duplicate entries
157 }
158 } // for each columm index in lclRow
159 } // if lclDiagCol is valid
160 } // numEnt != 0
161 } // for each lclRow
162
163 return diagCount;
164 } // if-else
165 }
166
168 template<class LO, class GO, class NT>
169 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
171 {
172 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
173
174 const auto rowMap = G.getRowMap ();
175 const auto colMap = G.getColMap ();
176 if (rowMap.get () == nullptr || colMap.get () == nullptr) {
177 return 0; // this process does not participate
178 }
179 else {
180 using inds_type = typename ::Tpetra::RowGraph<LO,GO,NT>::nonconst_local_inds_host_view_type;
181 inds_type lclColIndsBuf("lclColIndsBuf",G.getLocalMaxNumRowEntries());
182 const LO lclNumRows = static_cast<LO> (G.getLocalNumRows ());
183
184 LO diagCount = 0;
185 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
186 size_t numEntSizeT = G.getNumEntriesInLocalRow (lclRow);
187 const LO numEnt = static_cast<LO> (numEntSizeT);
188
189 inds_type lclColInds = Kokkos::subview(lclColIndsBuf,std::make_pair(0,numEnt));
190 G.getLocalRowCopy (lclRow, lclColInds, numEntSizeT);
191
192 if (numEnt != 0) {
193 const LO lclDiagCol =
194 colMap->getLocalElement (rowMap->getGlobalElement (lclRow));
195 // If it's not in the column Map, then there's no diagonal entry.
196 if (lclDiagCol != LOT::invalid ()) {
197 // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize
198 // for the sorted case.
199 for (LO k = 0; k < numEnt; ++k) {
200 if (lclDiagCol == lclColInds[k]) {
201 ++diagCount;
202 break; // don't count duplicate entries
203 }
204 } // for each columm index in lclRow
205 } // if lclDiagCol is valid
206 } // numEnt != 0
207 } // for each lclRow
208
209 return diagCount;
210 } // if-else
211 }
212
214 template<class LO, class GO, class NT>
215 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
217 {
218 const auto rowMap = G.getRowMap ();
219 if (rowMap.get () == nullptr) {
220 return 0; // this process does not participate
221 }
222 else {
223 typename ::Tpetra::RowGraph<LO,GO,NT>::global_inds_host_view_type
224 gblColInds;
225 const LO lclNumRows = static_cast<LO> (G.getLocalNumRows ());
226
227 LO diagCount = 0;
228 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
229 const GO gblRow = rowMap->getGlobalElement (lclRow);
230 G.getGlobalRowView (gblRow, gblColInds);
231 const LO numEnt = static_cast<LO> (gblColInds.size ());
232 if (numEnt != 0) {
233 // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize for
234 // the sorted case.
235 for (LO k = 0; k < numEnt; ++k) {
236 if (gblRow == gblColInds[k]) {
237 ++diagCount;
238 break; // don't count duplicate entries
239 }
240 } // for each column index in lclRow
241 } // if numEnt != 0
242 } // for each lclRow
243
244 return diagCount;
245 } // if-else
246 }
247
249 template<class LO, class GO, class NT>
250 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
252 {
253 using gids_type = typename ::Tpetra::RowGraph<LO,GO,NT>::nonconst_global_inds_host_view_type ;
254 const auto rowMap = G.getRowMap ();
255 if (rowMap.get () == nullptr) {
256 return 0; // this process does not participate
257 }
258 else {
259 gids_type gblColIndsBuf;
260 const LO lclNumRows = static_cast<LO> (G.getLocalNumRows ());
261
262 LO diagCount = 0;
263 for (LO lclRow = 0; lclRow < lclNumRows; ++lclRow) {
264 size_t numEntSizeT = G.getNumEntriesInLocalRow (lclRow);
265 const LO numEnt = static_cast<LO> (numEntSizeT);
266 if (static_cast<LO> (gblColIndsBuf.size ()) < numEnt) {
267 Kokkos::resize(gblColIndsBuf,numEnt);
268 }
269
270 gids_type gblColInds = Kokkos::subview(gblColIndsBuf,std::make_pair((LO)0, numEnt));
271 const GO gblRow = rowMap->getGlobalElement (lclRow);
272 G.getGlobalRowCopy (gblRow, gblColInds, numEntSizeT);
273
274 if (numEnt != 0) {
275 // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize for
276 // the sorted case.
277 for (LO k = 0; k < numEnt; ++k) {
278 if (gblRow == gblColInds[k]) {
279 ++diagCount;
280 break; // don't count duplicate entries
281 }
282 } // for each column index in lclRow
283 } // if numEnt != 0
284 } // for each lclRow
285
286 return diagCount;
287 } // if-else
288 }
289
296 template<class MatrixType>
298 static typename MatrixType::local_ordinal_type
299 getLocalNumDiags (const MatrixType& A)
300 {
301 using LO = typename MatrixType::local_ordinal_type;
302 using GO = typename MatrixType::global_ordinal_type;
303 using NT = typename MatrixType::node_type;
304 using row_graph_type = ::Tpetra::RowGraph<LO, GO, NT>;
305
306 auto G = A.getGraph ();
307 if (G.get () == nullptr) {
308 return 0;
309 }
310 else {
311 return GetLocalNumDiags<row_graph_type>::getLocalNumDiags (*G);
312 }
313 }
314 };
315
317 template<class LO, class GO, class NT>
318 struct GetLocalNumDiags< ::Tpetra::RowGraph<LO, GO, NT> > {
319 static LO
320 getLocalNumDiags (const ::Tpetra::RowGraph<LO, GO, NT>& G)
321 {
323
324 const crs_graph_type* G_crs = dynamic_cast<const crs_graph_type*> (&G);
325 if (G_crs != nullptr && G_crs->isFillComplete ()) {
326 return countLocalNumDiagsInFillCompleteGraph (*G_crs);
327 }
328 else {
329 if (G.isLocallyIndexed ()) {
330 if (G.supportsRowViews ()) {
332 }
333 else {
335 }
336 }
337 else if (G.isGloballyIndexed ()) {
338 if (G.supportsRowViews ()) {
340 }
341 else {
343 }
344 }
345 else { // G is empty
346 return 0;
347 }
348 }
349 }
350 };
351
353 template<class LO, class GO, class NT>
354 struct GetLocalNumDiags< ::Tpetra::CrsGraph<LO, GO, NT> > {
355 static LO
356 getLocalNumDiags (const ::Tpetra::CrsGraph<LO, GO, NT>& G)
357 {
358 using row_graph_type = ::Tpetra::RowGraph<LO, GO, NT>;
359 return GetLocalNumDiags<row_graph_type>::getLocalNumDiags (G);
360 }
361 };
362} // namespace Impl
363
366template<class CrsGraphType>
367typename CrsGraphType::local_ordinal_type
368getLocalNumDiags (const CrsGraphType& G)
369{
370 return Impl::GetLocalNumDiags<CrsGraphType>::getLocalNumDiags (G);
371}
372
375template<class CrsGraphType>
376typename CrsGraphType::global_ordinal_type
377getGlobalNumDiags (const CrsGraphType& G)
378{
379 using GO = typename CrsGraphType::global_ordinal_type;
380
381 const auto map = G.getRowMap ();
382 if (map.get () == nullptr) {
383 return GO (0); // this process should not participate
384 }
385 else {
386 const auto comm = map->getComm ();
387 if (comm.get () == nullptr) {
388 return GO (0); // this process should not participate
389 }
390 else {
391 const GO lclNumDiags = static_cast<GO> (getLocalNumDiags (G));
392
393 using Teuchos::REDUCE_SUM;
394 using Teuchos::reduceAll;
395 using Teuchos::outArg;
396
397 GO gblNumDiags {0};
398 reduceAll<int, GO> (*comm, REDUCE_SUM, lclNumDiags, outArg (gblNumDiags));
399 return gblNumDiags;
400 }
401 }
402}
403
404} // namespace Details
405} // namespace Tpetra
406
407#endif // TPETRA_DETAILS_GETNUMDIAGS_HPP
408
Import KokkosSparse::OrdinalTraits, a traits class for "invalid" (flag) values of integer types,...
typename::Tpetra::RowGraph< LO, GO, NT >::local_ordinal_type countLocalNumDiagsInNonFillCompleteGloballyIndexedGraphWithoutRowViews(const ::Tpetra::RowGraph< LO, GO, NT > &G)
Return local number of diagonal entries.
typename::Tpetra::RowGraph< LO, GO, NT >::local_ordinal_type countLocalNumDiagsInNonFillCompleteGloballyIndexedGraphWithRowViews(const ::Tpetra::RowGraph< LO, GO, NT > &G)
Return local number of diagonal entries.
MapType::local_ordinal_type getLocalDiagonalColumnIndex(const typename MapType::local_ordinal_type lclRow, const MapType &rowMap, const MapType &colMap)
Local columm index of diagonal entry.
typename::Tpetra::RowGraph< LO, GO, NT >::local_ordinal_type countLocalNumDiagsInNonFillCompleteLocallyIndexedGraphWithoutRowViews(const ::Tpetra::RowGraph< LO, GO, NT > &G)
Return local number of diagonal entries.
typename::Tpetra::RowGraph< LO, GO, NT >::local_ordinal_type countLocalNumDiagsInNonFillCompleteLocallyIndexedGraphWithRowViews(const ::Tpetra::RowGraph< LO, GO, NT > &G)
Return local number of diagonal entries.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
Teuchos::RCP< const map_type > getRowMap() const override
Returns the Map that describes the row distribution in this graph.
bool isFillComplete() const override
Whether the matrix is fill complete.
KOKKOS_INLINE_FUNCTION void operator()(const typename LocalMapType::local_ordinal_type lclRow, result_type &diagCount) const
Reduction function: result is (diagonal count, error count).
An abstract interface for graphs accessed by rows.
Nonmember function that computes a residual Computes R = B - A * X.
CrsGraphType::global_ordinal_type getGlobalNumDiags(const CrsGraphType &G)
Number of populated diagonal entries in the given sparse graph, over all processes in the graph's (MP...
CrsGraphType::local_ordinal_type getLocalNumDiags(const CrsGraphType &G)
Number of populated diagonal entries in the given sparse graph, on the calling (MPI) process.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Implementation of Tpetra::Details::getLocalNumDiags (see below).