Zoltan2
Loading...
Searching...
No Matches
Zoltan2_XpetraMultiVectorAdapter.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_XPETRAMULTIVECTORADAPTER_HPP_
15#define _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_
16
21
22#include <Xpetra_TpetraMultiVector.hpp>
23
24namespace Zoltan2 {
25
42template <typename User>
44public:
45
46#ifndef DOXYGEN_SHOULD_SKIP_THIS
48 typedef typename InputTraits<User>::impl_scalar_t impl_scalar_t;
49 typedef typename InputTraits<User>::lno_t lno_t;
50 typedef typename InputTraits<User>::gno_t gno_t;
51 typedef typename InputTraits<User>::part_t part_t;
52 typedef typename InputTraits<User>::node_t node_t;
53 typedef User user_t;
54 typedef User userCoord_t;
55
56 typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
57 typedef Xpetra::TpetraMultiVector<
58 scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
59#endif
60
76 XpetraMultiVectorAdapter(const RCP<const User> &invector,
77 std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides);
78
84 XpetraMultiVectorAdapter(const RCP<const User> &invector);
85
86
88 // The Adapter interface.
90
91 size_t getLocalNumIDs() const { return vector_->getLocalLength();}
92
93 void getIDsView(const gno_t *&ids) const
94 {
95 ids = map_->getLocalElementList().getRawPtr();
96 }
97
99 Kokkos::View<const gno_t *, typename node_t::device_type> &ids) const {
100 if (map_->lib() == Xpetra::UseTpetra) {
101 using device_type = typename node_t::device_type;
102 const xt_mvector_t *tvector =
103 dynamic_cast<const xt_mvector_t *>(vector_.get());
104 // MJ can be running Host, CudaSpace, or CudaUVMSpace while Map now
105 // internally never stores CudaUVMSpace so we may need a conversion.
106 // However Map stores both Host and CudaSpace so this could be improved
107 // if device_type was CudaSpace. Then we could add a new accessor to
108 // Map such as getMyGlobalIndicesDevice() which could be direct assigned
109 // here. Since Tpetra is still UVM dependent that is not going to happen
110 // yet so just leaving this as Host to device_type conversion for now.
111 ids = Kokkos::create_mirror_view_and_copy(device_type(),
112 tvector->getTpetra_MultiVector()->getMap()->getMyGlobalIndices());
113 }
114 else {
115 throw std::logic_error("getIDsKokkosView called but not on Tpetra!");
116 }
117 }
118
119 int getNumWeightsPerID() const { return numWeights_;}
120
121 void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
122 {
123 if(idx<0 || idx >= numWeights_)
124 {
125 std::ostringstream emsg;
126 emsg << __FILE__ << ":" << __LINE__
127 << " Invalid weight index " << idx << std::endl;
128 throw std::runtime_error(emsg.str());
129 }
130
131 size_t length;
132 weights_[idx].getStridedList(length, weights, stride);
133 }
134
135 void getWeightsKokkos2dView(Kokkos::View<scalar_t **,
136 typename node_t::device_type> &wgt) const {
137 typedef Kokkos::View<scalar_t**, typename node_t::device_type> view_t;
138 wgt = view_t("wgts", vector_->getLocalLength(), numWeights_);
139 typename view_t::host_mirror_type host_wgt = Kokkos::create_mirror_view(wgt);
140 for(int idx = 0; idx < numWeights_; ++idx) {
141 const scalar_t * weights;
142 size_t length;
143 int stride;
144 weights_[idx].getStridedList(length, weights, stride);
145 size_t fill_index = 0;
146 for(size_t n = 0; n < length; n += stride) {
147 host_wgt(fill_index++,idx) = weights[n];
148 }
149 }
150 Kokkos::deep_copy(wgt, host_wgt);
151 }
152
154 // The VectorAdapter interface.
156
157 int getNumEntriesPerID() const {return vector_->getNumVectors();}
158
159 void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const;
160
162 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
163 Kokkos::View<impl_scalar_t **, Kokkos::LayoutLeft,
164 typename node_t::device_type> & elements) const;
165
166 template <typename Adapter>
167 void applyPartitioningSolution(const User &in, User *&out,
168 const PartitioningSolution<Adapter> &solution) const;
169
170 template <typename Adapter>
171 void applyPartitioningSolution(const User &in, RCP<User> &out,
172 const PartitioningSolution<Adapter> &solution) const;
173
174private:
175
176 RCP<const User> invector_;
177 RCP<const x_mvector_t> vector_;
178 RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
179
180 int numWeights_;
181 ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
182};
183
185// Definitions
187
188template <typename User>
190 const RCP<const User> &invector,
191 std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
192 invector_(invector), vector_(), map_(),
193 numWeights_(weights.size()), weights_(weights.size())
194{
195 typedef StridedData<lno_t, scalar_t> input_t;
196
197 try {
198 RCP<x_mvector_t> tmp =
199 XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
200 vector_ = rcp_const_cast<const x_mvector_t>(tmp);
201 }
203
204 map_ = vector_->getMap();
205
206 size_t length = vector_->getLocalLength();
207
208 if (length > 0 && numWeights_ > 0){
209 int stride = 1;
210 for (int w=0; w < numWeights_; w++){
211 if (weightStrides.size())
212 stride = weightStrides[w];
213 ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*length, false);
214 weights_[w] = input_t(wgtV, stride);
215 }
216 }
217}
218
219
221template <typename User>
223 const RCP<const User> &invector):
224 invector_(invector), vector_(), map_(),
225 numWeights_(0), weights_()
226{
227 try {
228 RCP<x_mvector_t> tmp =
229 XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
230 vector_ = rcp_const_cast<const x_mvector_t>(tmp);
231 }
233
234 map_ = vector_->getMap();
235}
236
238template <typename User>
240 const scalar_t *&elements, int &stride, int idx) const
241{
242 size_t vecsize;
243 stride = 1;
244 elements = NULL;
245 if (map_->lib() == Xpetra::UseTpetra){
246 const xt_mvector_t *tvector =
247 dynamic_cast<const xt_mvector_t *>(vector_.get());
248
249 vecsize = tvector->getLocalLength();
250 if (vecsize > 0){
251 ArrayRCP<const scalar_t> data = tvector->getData(idx);
252 elements = data.get();
253 }
254 }
255 else{
256 throw std::logic_error("invalid underlying lib");
257 }
258}
259
261template <typename User>
263 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
264 Kokkos::View<impl_scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> & elements) const
265{
266 if (map_->lib() == Xpetra::UseTpetra){
267 const xt_mvector_t *tvector =
268 dynamic_cast<const xt_mvector_t *>(vector_.get());
269 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
270 auto view2d =
271 tvector->getTpetra_MultiVector()->template getLocalView<typename node_t::device_type>(Tpetra::Access::ReadWrite);
272 elements = view2d;
273 // CMS/KDD: Look at this stuff right here. Compare against a non-cuda build OR, look at core/driver/driverinputs/kuberry/kuberry.coords
274 // Ca try changing the kuberry.xml to use "input adapter" "BasicVector" rather than "XpetraMultiVector"
275
276 }
277 else {
278 throw std::logic_error("getEntriesKokkosView called but not using Tpetra!");
279 }
280}
281
283template <typename User>
284 template <typename Adapter>
286 const User &in, User *&out,
287 const PartitioningSolution<Adapter> &solution) const
288{
289 // Get an import list (rows to be received)
290 size_t numNewRows;
291 ArrayRCP<gno_t> importList;
292 try{
293 numNewRows = Zoltan2::getImportList<Adapter,
295 (solution, this, importList);
296 }
298
299 // Move the rows, creating a new vector.
300 RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewRows,
301 importList.getRawPtr());
302 out = outPtr.get();
303 outPtr.release();
304}
305
307template <typename User>
308 template <typename Adapter>
310 const User &in, RCP<User> &out,
311 const PartitioningSolution<Adapter> &solution) const
312{
313 // Get an import list (rows to be received)
314 size_t numNewRows;
315 ArrayRCP<gno_t> importList;
316 try{
317 numNewRows = Zoltan2::getImportList<Adapter,
319 (solution, this, importList);
320 }
322
323 // Move the rows, creating a new vector.
324 out = XpetraTraits<User>::doMigration(in, numNewRows,
325 importList.getRawPtr());
326}
327
328} //namespace Zoltan2
329
330#endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition Metric.cpp:39
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines the VectorAdapter interface.
Traits of Xpetra classes, including migration method.
typename BaseAdapter< User >::scalar_t scalar_t
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_t
typename InputTraits< User >::gno_t gno_t
typename InputTraits< User >::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
void getWeightsKokkos2dView(Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor.
void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
int getNumEntriesPerID() const
Return the number of vectors.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
void getEntriesKokkosView(Kokkos::View< impl_scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &elements) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
map_t::global_ordinal_type gno_t
Created by mbenlioglu on Aug 31, 2020.
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
static ArrayRCP< ArrayRCP< zscalar_t > > weights
The traits required of User input classes or structures.
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.
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.