Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_getGraphOffRankOffsets_decl.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_GETGRAPHOFFRANKOFFSETS_DECL_HPP
11#define TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DECL_HPP
12
17
18#include "TpetraCore_config.h"
19#include "Kokkos_Core.hpp"
20#include "Kokkos_StaticCrsGraph.hpp"
22#include <type_traits>
23
24namespace Tpetra {
25namespace Details {
26namespace Impl {
27
41template<class LO,
42 class GO,
43 class DeviceType,
44 class OffsetType = size_t>
46public:
47 typedef typename DeviceType::device_type device_type;
48 typedef OffsetType offset_type;
49 typedef ::Kokkos::View<offset_type*,
50 device_type,
51 ::Kokkos::MemoryUnmanaged> offsets_type;
52 typedef ::Kokkos::StaticCrsGraph<LO,
53 ::Kokkos::LayoutLeft,
54 device_type,
55 void, size_t> local_graph_type;
56 typedef ::Tpetra::Details::LocalMap<LO, GO, device_type> local_map_type;
57 typedef ::Kokkos::View<const typename local_graph_type::size_type*,
58 ::Kokkos::LayoutLeft,
59 device_type,
60 ::Kokkos::MemoryUnmanaged> row_offsets_type;
61 // This is unmanaged for performance, because we need to take
62 // subviews inside the functor.
63 typedef ::Kokkos::View<const LO*,
64 ::Kokkos::LayoutLeft,
65 device_type,
66 ::Kokkos::MemoryUnmanaged> lcl_col_inds_type;
67
69 GetGraphOffRankOffsets (const offsets_type& OffRankOffsets,
70 const local_map_type& lclColMap,
71 const local_map_type& lclDomMap,
72 const row_offsets_type& ptr,
73 const lcl_col_inds_type& ind);
74
76 KOKKOS_FUNCTION void operator() (const LO& lclRowInd) const;
77
78private:
79 offsets_type OffRankOffsets_;
80 local_map_type lclColMap_;
81 local_map_type lclDomMap_;
82 row_offsets_type ptr_;
83 lcl_col_inds_type ind_;
84 LO lclNumRows_;
85};
86
87} // namespace Impl
88
89template<class OffsetsType,
90 class LclMapType,
91 class LclGraphType>
92void
93getGraphOffRankOffsets (const OffsetsType& OffRankOffsets,
94 const LclMapType& lclColMap,
95 const LclMapType& lclDomMap,
96 const LclGraphType& lclGraph)
97{
98 typedef typename OffsetsType::non_const_value_type offset_type;
99 typedef typename LclMapType::local_ordinal_type LO;
100 typedef typename LclMapType::global_ordinal_type GO;
101 typedef typename LclMapType::device_type DT;
102
104
105 // The functor's constructor runs the functor.
106 impl_type impl (OffRankOffsets, lclColMap, lclDomMap, lclGraph.row_map, lclGraph.entries);
107}
108
109} // namespace Details
110} // namespace Tpetra
111
112#endif // TPETRA_DETAILS_GETGRAPHOFFRANKOFFSETS_DECL_HPP
Declaration and definition of the Tpetra::Map class, an implementation detail of Tpetra::Map.
Implementation detail of Tpetra::Details::getGraphOffRankOffsets, which in turn is an implementation ...
GetGraphOffRankOffsets(const offsets_type &OffRankOffsets, const local_map_type &lclColMap, const local_map_type &lclDomMap, const row_offsets_type &ptr, const lcl_col_inds_type &ind)
Constructor; also runs the functor.
KOKKOS_FUNCTION void operator()(const LO &lclRowInd) const
Kokkos::parallel_for loop body.
Nonmember function that computes a residual Computes R = B - A * X.
Namespace Tpetra contains the class and methods constituting the Tpetra library.