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
14#ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
15#define _ZOLTAN2_BASICVECTORADAPTER_HPP_
16
17#include "Zoltan2_Adapter.hpp"
20
21namespace Zoltan2 {
22
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
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
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
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 override { return numIds_; }
198
199 void getIDsView(const gno_t *&ids) const override { ids = idList_; }
200
201 void getIDsKokkosView(typename Base::ConstIdsDeviceView &ids) const override {
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 override { 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 override {
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 auto kWeights = this->kWeights_;
258 Kokkos::parallel_for(
259 size, KOKKOS_LAMBDA(const int id) {
260 deviceWgts(id) = kWeights(id, idx);
261 });
262
263 Kokkos::fence();
264 }
265
266 void
267 getWeightsDeviceView(typename Base::WeightsDeviceView &wgts) const override {
268 wgts = kWeights_;
269 }
270
272 // The VectorAdapter interface.
274
275 int getNumEntriesPerID() const override { return numEntriesPerID_; }
276
277 void getEntriesView(const scalar_t *&entries, int &stride,
278 int idx = 0) const override {
279 if (idx < 0 || idx >= numEntriesPerID_) {
280 std::ostringstream emsg;
281 emsg << __FILE__ << ":" << __LINE__ << " Invalid vector index " << idx
282 << std::endl;
283 throw std::runtime_error(emsg.str());
284 }
285 size_t length;
286 entries_[idx].getStridedList(length, entries, stride);
287 }
288
290 typename AdapterWithCoords<User>::CoordsDeviceView &entries) const override {
291 entries = kEntries_;
292 }
293
295 &entries) const override {
296 auto hostEntries = Kokkos::create_mirror_view(kEntries_);
297 Kokkos::deep_copy(hostEntries, kEntries_);
298 entries = hostEntries;
299 }
300
302 &entries) const override {
303 entries = kEntries_;
304 }
305
306private:
307 lno_t numIds_;
308 const gno_t *idList_;
309
310 int numEntriesPerID_;
311 int numWeights_;
312
313 // Old API variable members
314 ArrayRCP<StridedData<lno_t, scalar_t>> entries_;
315 ArrayRCP<StridedData<lno_t, scalar_t>> weights_;
316
317 // New API variable members
318 typename Base::IdsDeviceView kIds_;
320 typename Base::WeightsDeviceView kWeights_;
321
322 void createBasicVector(std::vector<const scalar_t *> &entries,
323 std::vector<int> &entryStride,
324 std::vector<const scalar_t *> &weights,
325 std::vector<int> &weightStrides) {
326 typedef StridedData<lno_t, scalar_t> input_t;
327
328 if (numIds_) {
329 // make kokkos ids
330 kIds_ = typename Base::IdsDeviceView(
331 Kokkos::ViewAllocateWithoutInitializing("ids"), numIds_);
332 auto host_kIds_ = Kokkos::create_mirror_view(kIds_);
333 for (int n = 0; n < numIds_; ++n) {
334 host_kIds_(n) = idList_[n];
335 }
336 Kokkos::deep_copy(kIds_, host_kIds_);
337
338 // make coordinates
339 int stride = 1;
340 entries_ = arcp(new input_t[numEntriesPerID_], 0, numEntriesPerID_, true);
341 for (int v = 0; v < numEntriesPerID_; v++) {
342 if (entryStride.size())
343 stride = entryStride[v];
344 ArrayRCP<const scalar_t> eltV(entries[v], 0, stride * numIds_, false);
345 entries_[v] = input_t(eltV, stride);
346 }
347
348 // setup kokkos entries
350 Kokkos::ViewAllocateWithoutInitializing("entries"), numIds_,
351 numEntriesPerID_);
352
353 size_t length;
354 const scalar_t *entriesPtr;
355
356 auto host_kokkos_entries = Kokkos::create_mirror_view(this->kEntries_);
357
358 for (int idx = 0; idx < numEntriesPerID_; idx++) {
359 entries_[idx].getStridedList(length, entriesPtr, stride);
360 size_t fill_index = 0;
361 for (int n = 0; n < numIds_; ++n) {
362 host_kokkos_entries(fill_index++, idx) = entriesPtr[n * stride];
363 }
364 }
365 Kokkos::deep_copy(this->kEntries_, host_kokkos_entries);
366 }
367
368 if (numWeights_) {
369 int stride = 1;
370 weights_ = arcp(new input_t[numWeights_], 0, numWeights_, true);
371 for (int w = 0; w < numWeights_; w++) {
372 if (weightStrides.size())
373 stride = weightStrides[w];
374 ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride * numIds_, false);
375 weights_[w] = input_t(wgtV, stride);
376 }
377
378 // set up final view with weights
379 kWeights_ = typename Base::WeightsDeviceView(
380 Kokkos::ViewAllocateWithoutInitializing("kokkos weights"), numIds_,
381 numWeights_);
382
383 // setup weights
384 auto host_weight_temp_values =
385 Kokkos::create_mirror_view(this->kWeights_);
386 for (int idx = 0; idx < numWeights_; ++idx) {
387 const scalar_t *weightsPtr;
388 size_t length;
389 weights_[idx].getStridedList(length, weightsPtr, stride);
390 size_t fill_index = 0;
391 for (size_t n = 0; n < length; n += stride) {
392 host_weight_temp_values(fill_index++, idx) = weightsPtr[n];
393 }
394 }
395 Kokkos::deep_copy(this->kWeights_, host_weight_temp_values);
396 }
397 }
398};
399
400} // namespace Zoltan2
401
402#endif
This file defines the StridedData class.
Defines the VectorAdapter interface.
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
typename BaseAdapter< User >::scalar_t scalar_t
typename CoordsDeviceView::host_mirror_type CoordsHostView
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_t
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
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 getIDsKokkosView(typename Base::ConstIdsDeviceView &ids) const override
void getWeightsHostView(typename Base::WeightsHostView &wgts) const override
void getWeightsView(const scalar_t *&weights, int &stride, int idx) 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 override
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
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const override
Provide a pointer to the elements of the specified vector.
int getNumWeightsPerID() const override
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.
int getNumEntriesPerID() const override
Return the number of vectors.
void getIDsView(const gno_t *&ids) const override
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.
virtual void getWeightsKokkos2dView(typename Base::WeightsDeviceView &wgt) const
size_t getLocalNumIDs() const override
Returns the number of objects on this process.
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...
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
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.