10#ifndef TPETRA_DETAILS_GETNUMDIAGS_HPP
11#define TPETRA_DETAILS_GETNUMDIAGS_HPP
20#include "Tpetra_CrsGraph.hpp"
21#include "Teuchos_CommHelpers.hpp"
33 template<
class LocalGraphType,
class LocalMapType>
39 G_ (
G), rowMap_ (rowMap), colMap_ (colMap)
44 using result_type =
typename LocalMapType::local_ordinal_type;
49 result_type& diagCount)
const
51 using LO =
typename LocalMapType::local_ordinal_type;
52 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
60 const LO
lclDiagCol = colMap_.getLocalElement (rowMap_.getGlobalElement (
lclRow));
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)
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;
89 using execution_space =
typename crs_graph_type::device_type::execution_space;
90 using policy_type = Kokkos::RangePolicy<execution_space, LO>;
92 const auto rowMap =
G.getRowMap ();
93 const auto colMap =
G.getColMap ();
94 if (rowMap.get () ==
nullptr || colMap.get () ==
nullptr) {
99 functor_type f (G.getLocalGraphDevice (), rowMap->getLocalMap (), colMap->getLocalMap ());
100 Kokkos::parallel_reduce (policy_type (0, G.getLocalNumRows ()), f, lclNumDiags);
114 template<
class MapType>
115 typename MapType::local_ordinal_type
120 return colMap.getLocalElement (rowMap.getGlobalElement (
lclRow));
124 template<
class LO,
class GO,
class NT>
125 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
128 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
130 const auto rowMap =
G.getRowMap ();
131 const auto colMap =
G.getColMap ();
132 if (rowMap.get () ==
nullptr || colMap.get () ==
nullptr) {
137 (!
G.supportsRowViews (), std::logic_error,
"Not implemented!");
139 typename ::Tpetra::RowGraph<LO, GO, NT>::local_inds_host_view_type
141 const LO
lclNumRows =
static_cast<LO
> (
G.getLocalNumRows ());
148 const LO
lclDiagCol = colMap->getLocalElement (rowMap->getGlobalElement (
lclRow));
168 template<
class LO,
class GO,
class NT>
169 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
172 using LOT = typename ::Tpetra::Details::OrdinalTraits<LO>;
174 const auto rowMap =
G.getRowMap ();
175 const auto colMap =
G.getColMap ();
176 if (rowMap.get () ==
nullptr || colMap.get () ==
nullptr) {
180 using inds_type = typename ::Tpetra::RowGraph<LO,GO,NT>::nonconst_local_inds_host_view_type;
182 const LO
lclNumRows =
static_cast<LO
> (
G.getLocalNumRows ());
194 colMap->getLocalElement (rowMap->getGlobalElement (
lclRow));
214 template<
class LO,
class GO,
class NT>
215 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
218 const auto rowMap =
G.getRowMap ();
219 if (rowMap.get () ==
nullptr) {
223 typename ::Tpetra::RowGraph<LO,GO,NT>::global_inds_host_view_type
225 const LO
lclNumRows =
static_cast<LO
> (
G.getLocalNumRows ());
249 template<
class LO,
class GO,
class NT>
250 typename ::Tpetra::RowGraph<LO, GO, NT>::local_ordinal_type
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) {
260 const LO
lclNumRows =
static_cast<LO
> (
G.getLocalNumRows ());
296 template<
class MatrixType>
298 static typename MatrixType::local_ordinal_type
301 using LO =
typename MatrixType::local_ordinal_type;
302 using GO =
typename MatrixType::global_ordinal_type;
303 using NT =
typename MatrixType::node_type;
306 auto G =
A.getGraph ();
307 if (
G.get () ==
nullptr) {
317 template<
class LO,
class GO,
class NT>
320 getLocalNumDiags (const ::Tpetra::RowGraph<LO, GO, NT>&
G)
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);
329 if (
G.isLocallyIndexed ()) {
330 if (
G.supportsRowViews ()) {
337 else if (
G.isGloballyIndexed ()) {
338 if (
G.supportsRowViews ()) {
353 template<
class LO,
class GO,
class NT>
356 getLocalNumDiags (const ::Tpetra::CrsGraph<LO, GO, NT>&
G)
366template<
class CrsGraphType>
367typename CrsGraphType::local_ordinal_type
375template<
class CrsGraphType>
376typename CrsGraphType::global_ordinal_type
379 using GO =
typename CrsGraphType::global_ordinal_type;
381 const auto map =
G.getRowMap ();
382 if (
map.get () ==
nullptr) {
386 const auto comm =
map->getComm ();
387 if (comm.get () ==
nullptr) {
393 using Teuchos::REDUCE_SUM;
394 using Teuchos::reduceAll;
395 using Teuchos::outArg;
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.
Struct that holds views of the contents of a CrsMatrix.
Kokkos::parallel_reduce functor for counting the local number of diagonal entries in a sparse graph.
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.
Implementation details of Tpetra.
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).