14#ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
15#define _ZOLTAN2_BASICVECTORADAPTER_HPP_
53#ifndef DOXYGEN_SHOULD_SKIP_THIS
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;
91 values.push_back(entries);
92 strides.push_back(entryStride);
94 weightValues.push_back(wgts);
95 weightStrides.push_back(wgtStride);
98 createBasicVector(values, strides, weightValues, weightStrides);
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()),
133 createBasicVector(entries, entryStride,
weights, weightStrides);
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;
173 strides.push_back(xStride);
177 strides.push_back(yStride);
181 strides.push_back(zStride);
187 weightValues.push_back(wgts);
188 weightStrides.push_back(wgtStride);
190 createBasicVector(values, strides, weightValues, weightStrides);
206 auto hostIds = Kokkos::create_mirror_view(this->kIds_);
207 Kokkos::deep_copy(hostIds, this->kIds_);
224 "Invalid weight index.");
227 weights_[idx].getStridedList(length,
weights, stride);
231 int idx = 0)
const override {
233 "Invalid weight index.");
236 typename Base::WeightsDeviceView1D(
"weights", kWeights_.extent(0));
239 hostWgts = Kokkos::create_mirror_view(weightsDevice);
240 Kokkos::deep_copy(hostWgts, weightsDevice);
244 auto hostWeights = Kokkos::create_mirror_view(kWeights_);
245 Kokkos::deep_copy(hostWeights, kWeights_);
250 int idx = 0)
const override {
252 "Invalid weight index.");
254 const auto size = kWeights_.extent(0);
255 deviceWgts =
typename Base::WeightsDeviceView1D(
"weights", size);
257 auto kWeights = this->kWeights_;
258 Kokkos::parallel_for(
259 size, KOKKOS_LAMBDA(
const int id) {
260 deviceWgts(
id) = kWeights(
id, idx);
278 int idx = 0)
const override {
279 if (idx < 0 || idx >= numEntriesPerID_) {
280 std::ostringstream emsg;
281 emsg << __FILE__ <<
":" << __LINE__ <<
" Invalid vector index " << idx
283 throw std::runtime_error(emsg.str());
286 entries_[idx].getStridedList(length, entries, stride);
295 &entries)
const override {
296 auto hostEntries = Kokkos::create_mirror_view(kEntries_);
297 Kokkos::deep_copy(hostEntries, kEntries_);
298 entries = hostEntries;
302 &entries)
const override {
308 const gno_t *idList_;
310 int numEntriesPerID_;
314 ArrayRCP<StridedData<lno_t, scalar_t>> entries_;
315 ArrayRCP<StridedData<lno_t, scalar_t>> weights_;
318 typename Base::IdsDeviceView kIds_;
320 typename Base::WeightsDeviceView kWeights_;
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) {
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];
336 Kokkos::deep_copy(kIds_, host_kIds_);
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);
350 Kokkos::ViewAllocateWithoutInitializing(
"entries"), numIds_,
356 auto host_kokkos_entries = Kokkos::create_mirror_view(this->kEntries_);
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];
365 Kokkos::deep_copy(this->kEntries_, host_kokkos_entries);
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);
379 kWeights_ =
typename Base::WeightsDeviceView(
380 Kokkos::ViewAllocateWithoutInitializing(
"kokkos weights"), numIds_,
384 auto host_weight_temp_values =
385 Kokkos::create_mirror_view(this->kWeights_);
386 for (
int idx = 0; idx < numWeights_; ++idx) {
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];
395 Kokkos::deep_copy(this->kWeights_, host_weight_temp_values);
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__)