Tpetra parallel linear algebra Version of the Day
Tpetra_Details_localRowOffsets_def.hpp
Go to the documentation of this file.
1// @HEADER
2// ***********************************************************************
3//
4// Tpetra: Templated Linear Algebra Services Package
5// Copyright (2008) Sandia Corporation
6//
7// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8// the U.S. Government retains certain rights in this software.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// ************************************************************************
38// @HEADER
39
40#ifndef TPETRA_DETAILS_LOCALROWOFFSETS_DEF_HPP
41#define TPETRA_DETAILS_LOCALROWOFFSETS_DEF_HPP
42
46
47#include "Tpetra_CrsGraph.hpp"
48#include "Tpetra_RowGraph.hpp"
51#include "Kokkos_Core.hpp"
52
53namespace Tpetra {
54namespace Details {
55namespace Impl {
56
57template <class LO, class GO, class NT>
58std::pair<typename LocalRowOffsetsResult<NT>::offsets_type, size_t>
59localRowCounts (const RowGraph<LO, GO, NT>& G)
60{
61 using result_type = LocalRowOffsetsResult<NT>;
62 using offsets_type = typename result_type::offsets_type;
63 using offset_type = typename result_type::offset_type;
64
65 const LO lclNumRows (G.getNodeNumRows ());
66 offsets_type entPerRow;
67 if (lclNumRows != 0) {
68 using Kokkos::view_alloc;
69 using Kokkos::WithoutInitializing;
70 entPerRow =
71 offsets_type (view_alloc ("entPerRow", WithoutInitializing),
72 lclNumRows);
73 }
74 using host = Kokkos::DefaultHostExecutionSpace;
75 auto entPerRow_h = Kokkos::create_mirror_view (host (), entPerRow);
76
77 // Don't trust G.getNodeMaxNumRowEntries() unless G is fillComplete.
78 // Even then, I would rather this method didn't exist (since it adds
79 // state and imposes overhead on fillComplete), and it's easy to
80 // compute ourselves here.
81 size_t maxNumEnt = 0;
82 for (LO i = 0; i < lclNumRows; ++i) {
83 const size_t lclNumEnt = G.getNumEntriesInLocalRow (i);
84 entPerRow_h[i] = offset_type (lclNumEnt);
85 maxNumEnt = maxNumEnt < lclNumEnt ? lclNumEnt : maxNumEnt;
86 }
87 Kokkos::deep_copy (entPerRow, entPerRow_h);
88 return {entPerRow, maxNumEnt};
89}
90
91template <class LO, class GO, class NT>
92LocalRowOffsetsResult<NT>
93localRowOffsetsFromRowGraph (const RowGraph<LO, GO, NT>& G)
94{
95 using result_type = LocalRowOffsetsResult<NT>;
96 using offsets_type = typename result_type::offsets_type;
97 using offset_type = typename result_type::offset_type;
98
99 offsets_type entPerRow;
100 size_t maxNumEnt = 0;
101 {
102 auto result = localRowCounts (G);
103 entPerRow = result.first;
104 maxNumEnt = result.second;
105 }
106
107 const LO lclNumRows (G.getNodeNumRows ());
108 offsets_type ptr;
109 offset_type nnz = 0;
110 if (lclNumRows != 0) {
111 using Kokkos::view_alloc;
112 using Kokkos::WithoutInitializing;
113 ptr = offsets_type (view_alloc ("ptr", WithoutInitializing),
114 lclNumRows + 1);
116 nnz = computeOffsetsFromCounts (ptr, entPerRow);
117 }
118 return {ptr, nnz, maxNumEnt};
119}
120
121template <class LO, class GO, class NT>
122LocalRowOffsetsResult<NT>
123localRowOffsetsFromFillCompleteCrsGraph (const CrsGraph<LO, GO, NT>& G)
124{
125 using Kokkos::view_alloc;
126 using Kokkos::WithoutInitializing;
127 using result_type = LocalRowOffsetsResult<NT>;
128 using offsets_type = typename result_type::offsets_type;
129 using offset_type = typename result_type::offset_type;
130
131 auto G_lcl = G.getLocalGraphDevice ();
132 offsets_type ptr (view_alloc ("ptr", WithoutInitializing),
133 G_lcl.row_map.extent (0));
134 Kokkos::deep_copy (ptr, G_lcl.row_map);
135
136 const offset_type nnz = G.getNodeNumEntries ();
137 const size_t maxNumEnt = G.getNodeMaxNumRowEntries ();
138 return {ptr, nnz, maxNumEnt};
139}
140
141} // namespace Impl
142
143template <class LO, class GO, class NT>
144LocalRowOffsetsResult<NT>
146{
147 if (G.isFillComplete ()) {
148 using crs_graph_type = CrsGraph<LO, GO, NT>;
149 const crs_graph_type* G_crs =
150 dynamic_cast<const crs_graph_type*> (&G);
151 if (G_crs != nullptr) {
152 return Impl::localRowOffsetsFromFillCompleteCrsGraph (*G_crs);
153 }
154 }
155 return Impl::localRowOffsetsFromRowGraph (G);
156}
157
158} // namespace Details
159} // namespace Tpetra
160
161//
162// Explicit instantiation macros
163//
164// Must be expanded from within the Tpetra namespace!
165//
166#define TPETRA_DETAILS_LOCALROWOFFSETS_INSTANT(LO, GO, NT) \
167namespace Details { \
168namespace Impl { \
169 \
170template std::pair<LocalRowOffsetsResult<NT>::offsets_type, size_t> \
171localRowCounts (const RowGraph<LO, GO, NT>& G); \
172 \
173template LocalRowOffsetsResult<NT> \
174localRowOffsetsFromRowGraph (const RowGraph<LO, GO, NT>& G); \
175 \
176template LocalRowOffsetsResult<NT> \
177localRowOffsetsFromFillCompleteCrsGraph (const CrsGraph<LO, GO, NT>& G); \
178 \
179} \
180 \
181template LocalRowOffsetsResult<NT> \
182localRowOffsets (const RowGraph<LO, GO, NT>& A); \
183}
184
185#endif // TPETRA_DETAILS_LOCALROWOFFSETS_DEF_HPP
Declare and define the functions Tpetra::Details::computeOffsetsFromCounts and Tpetra::computeOffsets...
Declaration and definition of Tpetra::Details::getEntryOnHost.
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
An abstract interface for graphs accessed by rows.
virtual bool isFillComplete() const =0
Whether fillComplete() has been called (without an intervening resumeFill()).
Implementation details of Tpetra.
LocalRowOffsetsResult< NT > localRowOffsets(const RowGraph< LO, GO, NT > &G)
Get local row offsets ("ptr", in compressed sparse row terms) for the given graph.
OffsetsViewType::non_const_value_type computeOffsetsFromCounts(const ExecutionSpace &execSpace, const OffsetsViewType &ptr, const CountsViewType &counts)
Compute offsets from counts.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.