Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Details_unpackCrsMatrixAndCombine_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_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
11#define TPETRA_DETAILS_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
12
13#include "TpetraCore_config.h"
15#include "Kokkos_DualView.hpp"
18#include "Tpetra_Details_DefaultTypes.hpp"
19
40
41#ifndef DOXYGEN_SHOULD_SKIP_THIS
42namespace Teuchos {
43// Forward declaration of Array
44template <class T>
45class Array;
46// Forward declaration of ArrayView
47template <class T>
48class ArrayView;
49} // namespace Teuchos
50#endif // DOXYGEN_SHOULD_SKIP_THIS
51
52namespace Tpetra {
53
54//
55// Users must never rely on anything in the Details namespace.
56//
57namespace Details {
58
97template <typename ST, typename LO, typename GO, typename NT>
99 const Teuchos::ArrayView<const char>& imports,
100 const Teuchos::ArrayView<const size_t>& numPacketsPerLID,
101 const Teuchos::ArrayView<const LO>& importLIDs,
102 size_t constantNumPackets,
104
105template <typename ST, typename LO, typename GO, typename NT>
106void unpackCrsMatrixAndCombineNew(
108 Kokkos::DualView<char*,
110 imports,
111 Kokkos::DualView<size_t*,
114 const Kokkos::DualView<const LO*,
116 const size_t constantNumPackets,
118
165//
174template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
175size_t
178 const Teuchos::ArrayView<const LocalOrdinal>& importLIDs,
179 const Teuchos::ArrayView<const char>& imports,
180 const Teuchos::ArrayView<const size_t>& numPacketsPerLID,
181 size_t constantNumPackets,
183 size_t numSameIDs,
184 const Teuchos::ArrayView<const LocalOrdinal>& permuteToLIDs,
185 const Teuchos::ArrayView<const LocalOrdinal>& permuteFromLIDs);
186
201
202template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
205 const Kokkos::View<LocalOrdinal const*,
206 Kokkos::Device<typename Node::device_type::execution_space,
207 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
209 ,
210 void, void
211#endif
212 >,
213 const Kokkos::View<const char*,
214 Kokkos::Device<typename Node::device_type::execution_space,
215 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
217 ,
218 void, void
219#endif
220 >,
221 const Kokkos::View<const size_t*,
222 Kokkos::Device<typename Node::device_type::execution_space,
223 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
225 ,
226 void, void
227#endif
228 >,
229 const size_t numSameIDs,
230 const Kokkos::View<LocalOrdinal const*,
231 Kokkos::Device<typename Node::device_type::execution_space,
232 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
234 ,
235 void, void
236#endif
237 >,
238 const Kokkos::View<LocalOrdinal const*,
239 Kokkos::Device<typename Node::device_type::execution_space,
240 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
242 ,
243 void, void
244#endif
245 >,
246 size_t TargetNumRows,
247 const int MyTargetPID,
248 Teuchos::ArrayRCP<size_t>& CRS_rowptr,
249 Teuchos::ArrayRCP<GlobalOrdinal>& CRS_colind,
250 Teuchos::ArrayRCP<Scalar>& CRS_vals,
251 const Teuchos::ArrayView<const int>& SourcePids,
252 Teuchos::Array<int>& TargetPids);
253
254template <typename Scalar, typename LocalOrdinal, typename GlobalOrdinal, typename Node>
257 const Kokkos::View<LocalOrdinal const*,
258 Kokkos::Device<typename Node::device_type::execution_space,
259 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
261 ,
262 void, void
263#endif
264 >,
265 const Kokkos::View<const char*,
266 Kokkos::Device<typename Node::device_type::execution_space,
267 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
269 ,
270 void, void
271#endif
272 >,
273 const Kokkos::View<const size_t*,
274 Kokkos::Device<typename Node::device_type::execution_space,
275 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
277 ,
278 void, void
279#endif
280 >,
281 const size_t numSameIDs,
282 const Kokkos::View<LocalOrdinal const*,
283 Kokkos::Device<typename Node::device_type::execution_space,
284 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
286 ,
287 void, void
288#endif
289 >,
290 const Kokkos::View<LocalOrdinal const*,
291 Kokkos::Device<typename Node::device_type::execution_space,
292 Tpetra::Details::DefaultTypes::comm_buffer_memory_space<typename Node::device_type>>
294 ,
295 void, void
296#endif
297 >,
298 size_t TargetNumRows,
299 const int MyTargetPID,
300 Kokkos::View<size_t*, typename Node::device_type>& /*crs_rowptr_d*/,
301 Kokkos::View<GlobalOrdinal*, typename Node::device_type>& /*crs_colind_d*/,
302 Kokkos::View<typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::impl_scalar_type*, typename Node::device_type>& /*crs_vals_d*/,
303 const Teuchos::ArrayView<const int>& SourcePids,
304 Kokkos::View<int*, typename Node::device_type>& /*TargetPids*/);
305
306} // namespace Details
307} // namespace Tpetra
308
309#endif // TPETRA_DETAILS_UNPACKCRSMATRIXANDCOMBINE_DECL_HPP
Declaration of Tpetra::CombineMode enum, and a function for setting a Tpetra::CombineMode parameter i...
Forward declaration of Tpetra::CrsMatrix.
Declaration of the Tpetra::DistObject class.
Struct that holds views of the contents of a CrsMatrix.
typename row_matrix_type::impl_scalar_type impl_scalar_type
The type used internally in place of Scalar.
Kokkos::Device< typename device_type::execution_space, buffer_memory_space > buffer_device_type
Kokkos::Device specialization for communication buffers.
Implementation details of Tpetra.
void unpackAndCombineIntoCrsArrays(const CrsGraph< LO, GO, NT > &sourceGraph, const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const typename CrsGraph< LO, GO, NT >::packet_type > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, const size_t constantNumPackets, const CombineMode combineMode, const size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs, size_t TargetNumRows, size_t TargetNumNonzeros, const int MyTargetPID, const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< GO > &CRS_colind, const Teuchos::ArrayView< const int > &SourcePids, Teuchos::Array< int > &TargetPids)
unpackAndCombineIntoCrsArrays
size_t unpackAndCombineWithOwningPIDsCount(const CrsGraph< LO, GO, NT > &sourceGraph, const Teuchos::ArrayView< const LO > &importLIDs, const Teuchos::ArrayView< const typename CrsGraph< LO, GO, NT >::packet_type > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, size_t constantNumPackets, CombineMode combineMode, size_t numSameIDs, const Teuchos::ArrayView< const LO > &permuteToLIDs, const Teuchos::ArrayView< const LO > &permuteFromLIDs)
Special version of Tpetra::Details::unpackCrsGraphAndCombine that also unpacks owning process ranks.
void unpackCrsMatrixAndCombine(const CrsMatrix< ST, LO, GO, NT > &sourceMatrix, const Teuchos::ArrayView< const char > &imports, const Teuchos::ArrayView< const size_t > &numPacketsPerLID, const Teuchos::ArrayView< const LO > &importLIDs, size_t constantNumPackets, CombineMode combineMode)
Unpack the imported column indices and values, and combine into matrix.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
CombineMode
Rule for combining data in an Import or Export.