Zoltan2
Loading...
Searching...
No Matches
Zoltan2_TpetraCrsGraphAdapter.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef _ZOLTAN2_TPETRACRSGRAPHADAPTER_HPP_
15#define _ZOLTAN2_TPETRACRSGRAPHADAPTER_HPP_
16
21#include <string>
22
23namespace Zoltan2 {
24
32template <typename User, typename UserCoord = User>
33class TpetraCrsGraphAdapter : public TpetraRowGraphAdapter<User, UserCoord> {
34
35public:
36#ifndef DOXYGEN_SHOULD_SKIP_THIS
39 using lno_t = typename InputTraits<User>::lno_t;
40 using gno_t = typename InputTraits<User>::gno_t;
41 using part_t = typename InputTraits<User>::part_t;
42 using node_t = typename InputTraits<User>::node_t;
43 using user_t = User;
44 using userCoord_t = UserCoord;
45
48#endif
49
59 TpetraCrsGraphAdapter(const RCP<const User> &graph, int nVtxWeights = 0,
60 int nEdgeWeights = 0);
61
62 void init(const RCP<const User> &graph);
63
66 RCP<const User> getUserGraph() const { return this->graph_; }
67
68 template <typename Adapter>
70 const User &in, User *&out,
71 const PartitioningSolution<Adapter> &solution) const;
72
73 template <typename Adapter>
75 const User &in, RCP<User> &out,
76 const PartitioningSolution<Adapter> &solution) const;
77};
78
80// Definitions
82
83template <typename User, typename UserCoord>
84typename TpetraCrsGraphAdapter<User, UserCoord>::Base::IdsDeviceView
85getColIds(const RCP<const User> &inmatrix) {
86 auto colIdsDevice = inmatrix->getLocalIndicesDevice();
87
88 auto colIdsGlobalDevice =
89 typename TpetraCrsGraphAdapter<User, UserCoord>::Base::IdsDeviceView("colIdsGlobalDevice", colIdsDevice.extent(0));
90 auto colMap = inmatrix->getColMap();
91 auto lclColMap = colMap->getLocalMap();
92
93 // Convert to global IDs using Tpetra::Map
94 Kokkos::parallel_for("colIdsGlobalDevice",
95 Kokkos::RangePolicy<typename User::node_type::execution_space>(
96 0, colIdsGlobalDevice.extent(0)),
97 KOKKOS_LAMBDA(const int i) {
98 colIdsGlobalDevice(i) =
99 lclColMap.getGlobalElement(colIdsDevice(i));
100 });
101
102 return colIdsGlobalDevice;
103}
104
105
106template <typename User, typename UserCoord>
107void TpetraCrsGraphAdapter<User, UserCoord>::init(const RCP<const User> &graph) {
108 auto colIdsDevice = graph->getLocalIndicesDevice();
109
110 auto colIdsGlobalDevice =
111 typename TpetraCrsGraphAdapter<User, UserCoord>::Base::IdsDeviceView("colIdsGlobalDevice", colIdsDevice.extent(0));
112 auto colMap = graph->getColMap();
113 auto lclColMap = colMap->getLocalMap();
114
115 // Convert to global IDs using Tpetra::Map
116 Kokkos::parallel_for("colIdsGlobalDevice",
117 Kokkos::RangePolicy<typename User::node_type::execution_space>(
118 0, colIdsGlobalDevice.extent(0)),
119 KOKKOS_LAMBDA(const int i) {
120 colIdsGlobalDevice(i) =
121 lclColMap.getGlobalElement(colIdsDevice(i));
122 });
123
124 this->adjIdsDevice_ = colIdsGlobalDevice;
125 this->offsDevice_ = graph->getLocalRowPtrsDevice();
126}
127
128template <typename User, typename UserCoord>
130 const RCP<const User> &graph, int nVtxWgts, int nEdgeWgts)
131 : TpetraRowGraphAdapter<User>(nVtxWgts, nEdgeWgts, graph) {
132
133 this->init(graph);
134
135 if (this->nWeightsPerVertex_ > 0) {
136
137 this->vertexWeightsDevice_ = typename Base::WeightsDeviceView(
138 "vertexWeightsDevice_", graph->getLocalNumRows(),
139 this->nWeightsPerVertex_);
140
141 this->vertexDegreeWeightsHost_ = typename Base::VtxDegreeHostView(
142 "vertexDegreeWeightsHost_", this->nWeightsPerVertex_);
143
144 for (int i = 0; i < this->nWeightsPerVertex_; ++i) {
145 this->vertexDegreeWeightsHost_(i) = false;
146 }
147 }
148
149 if (this->nWeightsPerEdge_) {
150 this->edgeWeightsDevice_ = typename Base::WeightsDeviceView(
151 "nWeightsPerEdge_", graph->getLocalNumRows(), this->nWeightsPerEdge_);
152 }
153}
154
156template <typename User, typename UserCoord>
157template <typename Adapter>
164
166template <typename User, typename UserCoord>
167template <typename Adapter>
174
175} // namespace Zoltan2
176
177#endif
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines TpetraRowGraphAdapter class.
Traits of Xpetra classes, including migration method.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< User >::gno_t gno_t
typename InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_t
GraphAdapter defines the interface for graph-based user data.
A PartitioningSolution is a solution to a partitioning problem.
Provides access for Zoltan2 to Tpetra::CrsGraph data.
RCP< const User > getUserGraph() const
Access to user's graph.
TpetraCrsGraphAdapter(const RCP< const User > &graph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
void init(const RCP< const User > &graph)
Provides access for Zoltan2 to Tpetra::RowGraph data.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Created by mbenlioglu on Aug 31, 2020.
TpetraCrsGraphAdapter< User, UserCoord >::Base::IdsDeviceView getColIds(const RCP< const User > &inmatrix)
default_offset_t offset_t
The data type to represent offsets.
default_gno_t gno_t
The ordinal type (e.g., int, long, int64_t) that can represent global counts and identifiers.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices.
default_part_t part_t
The data type to represent part numbers.
default_scalar_t scalar_t
The data type for weights and coordinates.