Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_TpetraCrsGraphAdapter.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 
50 #ifndef _ZOLTAN2_TPETRACRSGRAPHADAPTER_HPP_
51 #define _ZOLTAN2_TPETRACRSGRAPHADAPTER_HPP_
52 
54 #include <Zoltan2_StridedData.hpp>
56 #include <Zoltan2_XpetraTraits.hpp>
57 #include <string>
58 
59 namespace Zoltan2 {
60 
68 template <typename User, typename UserCoord = User>
69 class TpetraCrsGraphAdapter : public TpetraRowGraphAdapter<User, UserCoord> {
70 
71 public:
72 #ifndef DOXYGEN_SHOULD_SKIP_THIS
73  using scalar_t = typename InputTraits<User>::scalar_t;
74  using offset_t = typename InputTraits<User>::offset_t;
75  using lno_t = typename InputTraits<User>::lno_t;
76  using gno_t = typename InputTraits<User>::gno_t;
77  using part_t = typename InputTraits<User>::part_t;
78  using node_t = typename InputTraits<User>::node_t;
79  using user_t = User;
80  using userCoord_t = UserCoord;
81 
82  using Base = GraphAdapter<User, UserCoord>;
84 #endif
85 
95  TpetraCrsGraphAdapter(const RCP<const User> &graph, int nVtxWeights = 0,
96  int nEdgeWeights = 0);
97 
100  RCP<const User> getUserGraph() const { return this->graph_; }
101 
102  template <typename Adapter>
104  const User &in, User *&out,
105  const PartitioningSolution<Adapter> &solution) const;
106 
107  template <typename Adapter>
109  const User &in, RCP<User> &out,
110  const PartitioningSolution<Adapter> &solution) const;
111 };
112 
114 // Definitions
116 
117 template <typename User, typename UserCoord>
119  const RCP<const User> &graph, int nVtxWgts, int nEdgeWgts)
120  : TpetraRowGraphAdapter<User>(nVtxWgts, nEdgeWgts, graph) {
121  auto adjIdsHost = graph->getLocalIndicesHost();
122 
123  auto adjIdsGlobalHost =
124  typename Base::IdsHostView("adjIdsGlobalHost", adjIdsHost.extent(0));
125  auto colMap = graph->getColMap();
126 
127  // Convert to global IDs using Tpetra::Map
128  Kokkos::parallel_for("adjIdsGlobalHost",
129  Kokkos::RangePolicy<Kokkos::HostSpace::execution_space>(
130  0, adjIdsGlobalHost.extent(0)),
131  [=](const int i) {
132  adjIdsGlobalHost(i) =
133  colMap->getGlobalElement(adjIdsHost(i));
134  });
135 
136  auto adjIdsDevice = Kokkos::create_mirror_view_and_copy(
137  typename Base::device_t(), adjIdsGlobalHost);
138 
139  this->adjIdsDevice_ = adjIdsDevice;
140  this->offsDevice_ = graph->getLocalRowPtrsDevice();
141 
142  if (this->nWeightsPerVertex_ > 0) {
143 
144  this->vertexWeightsDevice_ = typename Base::WeightsDeviceView(
145  "vertexWeightsDevice_", graph->getLocalNumRows(),
146  this->nWeightsPerVertex_);
147 
148  this->vertexDegreeWeightsHost_ = typename Base::VtxDegreeHostView(
149  "vertexDegreeWeightsHost_", this->nWeightsPerVertex_);
150 
151  for (int i = 0; i < this->nWeightsPerVertex_; ++i) {
152  this->vertexDegreeWeightsHost_(i) = false;
153  }
154  }
155 
156  if (this->nWeightsPerEdge_) {
157  this->edgeWeightsDevice_ = typename Base::WeightsDeviceView(
158  "nWeightsPerEdge_", graph->getLocalNumRows(), this->nWeightsPerEdge_);
159  }
160 }
161 
163 template <typename User, typename UserCoord>
164 template <typename Adapter>
166  const User &in, User *&out,
167  const PartitioningSolution<Adapter> &solution) const {
169  solution);
170 }
171 
173 template <typename User, typename UserCoord>
174 template <typename Adapter>
176  const User &in, RCP<User> &out,
177  const PartitioningSolution<Adapter> &solution) const {
179  solution);
180 }
181 
182 } // namespace Zoltan2
183 
184 #endif
Helper functions for Partitioning Problems.
typename InputTraits< User >::scalar_t scalar_t
GraphAdapter defines the interface for graph-based user data.
Defines TpetraRowGraphAdapter class.
typename node_t::device_type device_t
default_part_t part_t
The data type to represent part numbers.
default_offset_t offset_t
The data type to represent offsets.
Provides access for Zoltan2 to Tpetra::CrsGraph data.
Traits of Xpetra classes, including migration method.
typename InputTraits< User >::part_t part_t
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
typename InputTraits< User >::node_t node_t
TpetraCrsGraphAdapter(const RCP< const User > &graph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
A PartitioningSolution is a solution to a partitioning problem.
typename InputTraits< User >::gno_t gno_t
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
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.
typename InputTraits< User >::offset_t offset_t
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
RCP< const User > getUserGraph() const
Access to user&#39;s graph.
typename InputTraits< User >::lno_t lno_t
default_scalar_t scalar_t
The data type for weights and coordinates.
Provides access for Zoltan2 to Tpetra::RowGraph data.
This file defines the StridedData class.