Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Tpetra_Details_determineLocalTriangularStructure.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_DETERMINELOCALTRIANGULARSTRUCTURE_HPP
11 #define TPETRA_DETAILS_DETERMINELOCALTRIANGULARSTRUCTURE_HPP
12 
19 
20 #include "Kokkos_Core.hpp"
22 
23 namespace Tpetra {
24 namespace Details {
25 
29 template<class LO>
39 };
40 
41 namespace Impl {
55  template<class LocalGraphType, class LocalMapType>
57  public:
58  // Result can't be more than the number of local rows, so
59  // local_ordinal_type is appropriate.
60  using result_type =
62 
72  DetermineLocalTriangularStructure (const LocalGraphType& G,
73  const LocalMapType& rowMap,
74  const LocalMapType& colMap,
75  const bool ignoreMapsForTriangularStructure) :
76  G_ (G),
77  rowMap_ (rowMap),
78  colMap_ (colMap),
79  ignoreMapsForTriangularStructure_ (ignoreMapsForTriangularStructure)
80  {}
81 
83  KOKKOS_INLINE_FUNCTION void init (result_type& dst) const
84  {
85  dst.diagCount = 0;
86  dst.maxNumRowEnt = 0;
87  dst.couldBeLowerTriangular = true; // well, we don't know yet, do we?
88  dst.couldBeUpperTriangular = true; // ditto
89  }
90 
91  KOKKOS_INLINE_FUNCTION void
92  join (result_type& dst,
93  const result_type& src) const
94  {
95  dst.diagCount += src.diagCount;
96  dst.maxNumRowEnt = (src.maxNumRowEnt > dst.maxNumRowEnt) ?
97  src.maxNumRowEnt : dst.maxNumRowEnt;
98  dst.couldBeLowerTriangular &= src.couldBeLowerTriangular;
99  dst.couldBeUpperTriangular &= src.couldBeUpperTriangular;
100  }
101 
103  KOKKOS_INLINE_FUNCTION void
104  operator () (const typename LocalMapType::local_ordinal_type lclRow,
105  result_type& result) const
106  {
107  using LO = typename LocalMapType::local_ordinal_type;
108  using GO = typename LocalMapType::global_ordinal_type;
109  using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
110 
111  auto G_row = G_.rowConst (lclRow);
112  const LO numEnt = G_row.length;
113  if (numEnt != 0) {
114  result.maxNumRowEnt = (numEnt > result.maxNumRowEnt) ?
115  numEnt : result.maxNumRowEnt;
116  // Use global row and column indices to find the diagonal
117  // entry. Caller promises that local row index is in the row
118  // Map on the calling process.
119  const GO gblDiagCol = rowMap_.getGlobalElement (lclRow);
120  const LO lclDiagCol = colMap_.getLocalElement (gblDiagCol);
121  // If it's not in the column Map, then there's no diagonal entry.
122  if (lclDiagCol != LOT::invalid ()) {
123  // TODO (mfh 25 Apr 2018) Use findRelOffset to optimize for
124  // the sorted case, but note that it requires operator[].
125  bool foundDiag = false; // don't count duplicates
126 
127  if (ignoreMapsForTriangularStructure_) {
128  for (LO k = 0; k < numEnt && ! foundDiag; ++k) {
129  const LO lclCol = G_row(k);
130  if (lclCol == lclDiagCol) {
131  foundDiag = true;
132  }
133  }
134  // mfh 30 Apr 2018: See GitHub Issue #2658. Per
135  // current Tpetra::CrsGraph::computeLocalConstants
136  // behavior, assume that local column indices are
137  // sorted in each row.
138  if (numEnt > LO (0)) {
139  const LO smallestLclCol = G_row(0);
140  const LO largestLclCol = G_row(numEnt-1); // could be same
141 
142  if (smallestLclCol < lclRow) {
143  result.couldBeUpperTriangular = false;
144  }
145  if (lclRow < largestLclCol) {
146  result.couldBeLowerTriangular = false;
147  }
148  }
149  }
150  else {
151  for (LO k = 0; k < numEnt &&
152  ((! foundDiag) ||
153  result.couldBeLowerTriangular ||
154  result.couldBeUpperTriangular);
155  ++k) {
156  const LO lclCol = G_row(k);
157  if (lclCol == lclDiagCol) {
158  foundDiag = true;
159  }
160  else {
161  const GO gblCol = colMap_.getGlobalElement (lclCol);
162  if (gblCol < gblDiagCol) {
163  result.couldBeUpperTriangular = false;
164  }
165  if (gblDiagCol < gblCol) {
166  result.couldBeLowerTriangular = false;
167  }
168  }
169  } // for each entry in lclRow
170  } // if-else ignoreMapsForTriangularStructure
171 
172  if (foundDiag) {
173  ++(result.diagCount);
174  }
175  }
176  }
177  }
178 
179  private:
180  LocalGraphType G_;
181  LocalMapType rowMap_;
182  LocalMapType colMap_;
183  bool ignoreMapsForTriangularStructure_;
184  };
185 
186 } // namespace Impl
187 
205 template<class LocalGraphType, class LocalMapType>
206 LocalTriangularStructureResult<typename LocalMapType::local_ordinal_type>
207 determineLocalTriangularStructure (const LocalGraphType& G,
208  const LocalMapType& rowMap,
209  const LocalMapType& colMap,
210  const bool ignoreMapsForTriangularStructure)
211 {
212  using LO = typename LocalMapType::local_ordinal_type;
213  using execution_space = typename LocalGraphType::device_type::execution_space;
214  using range_type = Kokkos::RangePolicy<execution_space, LO>;
215  using functor_type =
217 
218  LocalTriangularStructureResult<LO> result {0, 0, true, true};
219  Kokkos::parallel_reduce ("Tpetra::Details::determineLocalTriangularStructure",
220  range_type (0, G.numRows ()),
221  functor_type (G, rowMap, colMap,
222  ignoreMapsForTriangularStructure),
223  result);
224  return result;
225 }
226 
227 } // namespace Details
228 } // namespace Tpetra
229 
230 #endif // TPETRA_DETAILS_DETERMINELOCALTRIANGULARSTRUCTURE_HPP
Import KokkosSparse::OrdinalTraits, a traits class for &quot;invalid&quot; (flag) values of integer types...
bool couldBeLowerTriangular
Whether the graph is locally structurally lower triangular.
KOKKOS_INLINE_FUNCTION void init(result_type &dst) const
Set the initial value of the reduction result.
KOKKOS_INLINE_FUNCTION void operator()(const typename LocalMapType::local_ordinal_type lclRow, result_type &result) const
Reduction operator: result is (diagonal count, error count).
DetermineLocalTriangularStructure(const LocalGraphType &G, const LocalMapType &rowMap, const LocalMapType &colMap, const bool ignoreMapsForTriangularStructure)
Constructor.
LocalTriangularStructureResult< typename LocalMapType::local_ordinal_type > determineLocalTriangularStructure(const LocalGraphType &G, const LocalMapType &rowMap, const LocalMapType &colMap, const bool ignoreMapsForTriangularStructure)
Count the local number of diagonal entries in a local sparse graph, and determine whether the local p...
bool couldBeUpperTriangular
Whether the graph is locally structurally upper triangular.
Implementation of Tpetra::Details::determineLocalTriangularStructure (which see below).