Zoltan2
Loading...
Searching...
No Matches
Zoltan2_BasicVectorAdapter.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
13
14#ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
15#define _ZOLTAN2_BASICVECTORADAPTER_HPP_
16
17#include "Zoltan2_Adapter.hpp"
20
21namespace Zoltan2 {
22
49
50template <typename User> class BasicVectorAdapter : public VectorAdapter<User> {
51
52public:
53#ifndef DOXYGEN_SHOULD_SKIP_THIS
54
56 using lno_t = typename InputTraits<User>::lno_t;
57 using gno_t = typename InputTraits<User>::gno_t;
59 using part_t = typename InputTraits<User>::part_t;
60 using node_t = typename InputTraits<User>::node_t;
61 using user_t = User;
62 using Base = BaseAdapter<User>;
63
64#endif
65
80
81 BasicVectorAdapter(lno_t numIds, const gno_t *ids, const scalar_t *entries,
82 int entryStride = 1, bool usewgts = false,
83 const scalar_t *wgts = NULL, int wgtStride = 1)
84 : numIds_(numIds), idList_(ids), numEntriesPerID_(1),
85 numWeights_(usewgts == true) {
86 std::vector<const scalar_t *> values;
87 std::vector<int> strides;
88 std::vector<const scalar_t *> weightValues;
89 std::vector<int> weightStrides;
90
91 values.push_back(entries);
92 strides.push_back(entryStride);
93 if (usewgts) {
94 weightValues.push_back(wgts);
95 weightStrides.push_back(wgtStride);
96 }
97
98 createBasicVector(values, strides, weightValues, weightStrides);
99 }
100
125
126 BasicVectorAdapter(lno_t numIds, const gno_t *ids,
127 std::vector<const scalar_t *> &entries,
128 std::vector<int> &entryStride,
129 std::vector<const scalar_t *> &weights,
130 std::vector<int> &weightStrides)
131 : numIds_(numIds), idList_(ids), numEntriesPerID_(entries.size()),
132 numWeights_(weights.size()) {
133 createBasicVector(entries, entryStride, weights, weightStrides);
134 }
135
161
162 BasicVectorAdapter(lno_t numIds, const gno_t *ids, const scalar_t *x,
163 const scalar_t *y, const scalar_t *z, int xStride = 1,
164 int yStride = 1, int zStride = 1, bool usewgts = false,
165 const scalar_t *wgts = NULL, int wgtStride = 1)
166 : numIds_(numIds), idList_(ids), numEntriesPerID_(0),
167 numWeights_(usewgts == true) {
168 std::vector<const scalar_t *> values, weightValues;
169 std::vector<int> strides, weightStrides;
170
171 if (x) {
172 values.push_back(x);
173 strides.push_back(xStride);
174 numEntriesPerID_++;
175 if (y) {
176 values.push_back(y);
177 strides.push_back(yStride);
178 numEntriesPerID_++;
179 if (z) {
180 values.push_back(z);
181 strides.push_back(zStride);
182 numEntriesPerID_++;
183 }
184 }
185 }
186 if (usewgts) {
187 weightValues.push_back(wgts);
188 weightStrides.push_back(wgtStride);
189 }
190 createBasicVector(values, strides, weightValues, weightStrides);
191 }
192
194 // The Adapter interface.
196
197 size_t getLocalNumIDs() const { return numIds_; }
198
199 void getIDsView(const gno_t *&ids) const { ids = idList_; }
200
201 void getIDsKokkosView(typename Base::ConstIdsDeviceView &ids) const {
202 ids = this->kIds_;
203 }
204
205 void getIDsHostView(typename Base::ConstIdsHostView &ids) const override {
206 auto hostIds = Kokkos::create_mirror_view(this->kIds_);
207 Kokkos::deep_copy(hostIds, this->kIds_);
208 ids = hostIds;
209 }
210
211 void getIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override {
212 ids = this->kIds_;
213 }
214
215 int getNumWeightsPerID() const { return numWeights_; }
216
217 virtual void
218 getWeightsKokkos2dView(typename Base::WeightsDeviceView &wgt) const {
219 wgt = kWeights_;
220 }
221
222 void getWeightsView(const scalar_t *&weights, int &stride, int idx) const {
223 AssertCondition((idx >= 0) and (idx < numWeights_),
224 "Invalid weight index.");
225
226 size_t length;
227 weights_[idx].getStridedList(length, weights, stride);
228 }
229
230 void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts,
231 int idx = 0) const override {
232 AssertCondition((idx >= 0) and (idx < numWeights_),
233 "Invalid weight index.");
234
235 auto weightsDevice =
236 typename Base::WeightsDeviceView1D("weights", kWeights_.extent(0));
237 getWeightsDeviceView(weightsDevice, idx);
238
239 hostWgts = Kokkos::create_mirror_view(weightsDevice);
240 Kokkos::deep_copy(hostWgts, weightsDevice);
241 }
242
243 void getWeightsHostView(typename Base::WeightsHostView &wgts) const override {
244 auto hostWeights = Kokkos::create_mirror_view(kWeights_);
245 Kokkos::deep_copy(hostWeights, kWeights_);
246 wgts = hostWeights;
247 }
248
249 void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts,
250 int idx = 0) const override {
251 AssertCondition((idx >= 0) and (idx < numWeights_),
252 "Invalid weight index.");
253
254 const auto size = kWeights_.extent(0);
255 deviceWgts = typename Base::WeightsDeviceView1D("weights", size);
256
257 Kokkos::parallel_for(
258 size, KOKKOS_CLASS_LAMBDA(const int id) {
259 deviceWgts(id) = kWeights_(id, idx);
260 });
261
262 Kokkos::fence();
263 }
264
265 void
266 getWeightsDeviceView(typename Base::WeightsDeviceView &wgts) const override {
267 wgts = kWeights_;
268 }
269
271 // The VectorAdapter interface.
273
274 int getNumEntriesPerID() const { return numEntriesPerID_; }
275
276 void getEntriesView(const scalar_t *&entries, int &stride,
277 int idx = 0) const {
278 if (idx < 0 || idx >= numEntriesPerID_) {
279 std::ostringstream emsg;
280 emsg << __FILE__ << ":" << __LINE__ << " Invalid vector index " << idx
281 << std::endl;
282 throw std::runtime_error(emsg.str());
283 }
284 size_t length;
285 entries_[idx].getStridedList(length, entries, stride);
286 }
287
289 typename AdapterWithCoords<User>::CoordsDeviceView &entries) const {
290 entries = kEntries_;
291 }
292
294 &entries) const override {
295 auto hostEntries = Kokkos::create_mirror_view(kEntries_);
296 Kokkos::deep_copy(hostEntries, kEntries_);
297 entries = hostEntries;
298 }
299
301 &entries) const override {
302 entries = kEntries_;
303 }
304
305private:
306 lno_t numIds_;
307 const gno_t *idList_;
308
309 int numEntriesPerID_;
310 int numWeights_;
311
312 // Old API variable members
313 ArrayRCP<StridedData<lno_t, scalar_t>> entries_;
314 ArrayRCP<StridedData<lno_t, scalar_t>> weights_;
315
316 // New API variable members
317 typename Base::IdsDeviceView kIds_;
319 typename Base::WeightsDeviceView kWeights_;
320
321 void createBasicVector(std::vector<const scalar_t *> &entries,
322 std::vector<int> &entryStride,
323 std::vector<const scalar_t *> &weights,
324 std::vector<int> &weightStrides) {
325 typedef StridedData<lno_t, scalar_t> input_t;
326
327 if (numIds_) {
328 // make kokkos ids
329 kIds_ = typename Base::IdsDeviceView(
330 Kokkos::ViewAllocateWithoutInitializing("ids"), numIds_);
331 auto host_kIds_ = Kokkos::create_mirror_view(kIds_);
332 for (int n = 0; n < numIds_; ++n) {
333 host_kIds_(n) = idList_[n];
334 }
335 Kokkos::deep_copy(kIds_, host_kIds_);
336
337 // make coordinates
338 int stride = 1;
339 entries_ = arcp(new input_t[numEntriesPerID_], 0, numEntriesPerID_, true);
340 for (int v = 0; v < numEntriesPerID_; v++) {
341 if (entryStride.size())
342 stride = entryStride[v];
343 ArrayRCP<const scalar_t> eltV(entries[v], 0, stride * numIds_, false);
344 entries_[v] = input_t(eltV, stride);
345 }
346
347 // setup kokkos entries
349 Kokkos::ViewAllocateWithoutInitializing("entries"), numIds_,
350 numEntriesPerID_);
351
352 size_t length;
353 const scalar_t *entriesPtr;
354
355 auto host_kokkos_entries = Kokkos::create_mirror_view(this->kEntries_);
356
357 for (int idx = 0; idx < numEntriesPerID_; idx++) {
358 entries_[idx].getStridedList(length, entriesPtr, stride);
359 size_t fill_index = 0;
360 for (int n = 0; n < numIds_; ++n) {
361 host_kokkos_entries(fill_index++, idx) = entriesPtr[n * stride];
362 }
363 }
364 Kokkos::deep_copy(this->kEntries_, host_kokkos_entries);
365 }
366
367 if (numWeights_) {
368 int stride = 1;
369 weights_ = arcp(new input_t[numWeights_], 0, numWeights_, true);
370 for (int w = 0; w < numWeights_; w++) {
371 if (weightStrides.size())
372 stride = weightStrides[w];
373 ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride * numIds_, false);
374 weights_[w] = input_t(wgtV, stride);
375 }
376
377 // set up final view with weights
378 kWeights_ = typename Base::WeightsDeviceView(
379 Kokkos::ViewAllocateWithoutInitializing("kokkos weights"), numIds_,
380 numWeights_);
381
382 // setup weights
383 auto host_weight_temp_values =
384 Kokkos::create_mirror_view(this->kWeights_);
385 for (int idx = 0; idx < numWeights_; ++idx) {
386 const scalar_t *weightsPtr;
387 size_t length;
388 weights_[idx].getStridedList(length, weightsPtr, stride);
389 size_t fill_index = 0;
390 for (size_t n = 0; n < length; n += stride) {
391 host_weight_temp_values(fill_index++, idx) = weightsPtr[n];
392 }
393 }
394 Kokkos::deep_copy(this->kWeights_, host_weight_temp_values);
395 }
396 }
397};
398
399} // namespace Zoltan2
400
401#endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition Metric.cpp:39
This file defines the StridedData class.
Defines the VectorAdapter interface.
typename CoordsDeviceView::HostMirror CoordsHostView
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
typename BaseAdapter< User >::scalar_t scalar_t
typename InputTraits< User >::scalar_t scalar_t
BasicVectorAdapter(lno_t numIds, const gno_t *ids, std::vector< const scalar_t * > &entries, std::vector< int > &entryStride, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor for multivector (a set of vectors sharing the same global numbering and data distribution...
void getWeightsHostView(typename Base::WeightsHostView &wgts) const override
void getEntriesHostView(typename AdapterWithCoords< User >::CoordsHostView &entries) const override
Provide a Kokkos view (Host side) to the elements of the specified vector.
void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
void getEntriesKokkosView(typename AdapterWithCoords< User >::CoordsDeviceView &entries) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
void getWeightsDeviceView(typename Base::WeightsDeviceView &wgts) const override
void getIDsHostView(typename Base::ConstIdsHostView &ids) const override
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getEntriesDeviceView(typename AdapterWithCoords< User >::CoordsDeviceView &entries) const override
Provide a Kokkos view (Device side) to the elements of the specified vector.
BasicVectorAdapter(lno_t numIds, const gno_t *ids, const scalar_t *entries, int entryStride=1, bool usewgts=false, const scalar_t *wgts=NULL, int wgtStride=1)
Constructor for one vector with (optionally) one weight.
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
virtual void getWeightsKokkos2dView(typename Base::WeightsDeviceView &wgt) const
void getIDsKokkosView(typename Base::ConstIdsDeviceView &ids) const
void getIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override
BasicVectorAdapter(lno_t numIds, const gno_t *ids, const scalar_t *x, const scalar_t *y, const scalar_t *z, int xStride=1, int yStride=1, int zStride=1, bool usewgts=false, const scalar_t *wgts=NULL, int wgtStride=1)
A simple constructor for coordinate-based problems with dimension 1, 2 or 3 and (optionally) one weig...
void getIDsView(const gno_t *&ids) const
int getNumEntriesPerID() const
Return the number of vectors.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
Created by mbenlioglu on Aug 31, 2020.
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)
static ArrayRCP< ArrayRCP< zscalar_t > > weights
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.