50 #ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
51 #define _ZOLTAN2_BASICVECTORADAPTER_HPP_
89 #ifndef DOXYGEN_SHOULD_SKIP_THIS
118 int entryStride = 1,
bool usewgts =
false,
119 const scalar_t *wgts = NULL,
int wgtStride = 1)
120 : numIds_(numIds), idList_(ids), numEntriesPerID_(1),
121 numWeights_(usewgts == true) {
122 std::vector<const scalar_t *> values;
123 std::vector<int> strides;
124 std::vector<const scalar_t *> weightValues;
125 std::vector<int> weightStrides;
127 values.push_back(entries);
128 strides.push_back(entryStride);
130 weightValues.push_back(wgts);
131 weightStrides.push_back(wgtStride);
134 createBasicVector(values, strides, weightValues, weightStrides);
163 std::vector<const scalar_t *> &entries,
164 std::vector<int> &entryStride,
165 std::vector<const scalar_t *> &
weights,
166 std::vector<int> &weightStrides)
167 : numIds_(numIds), idList_(ids), numEntriesPerID_(entries.size()),
168 numWeights_(weights.size()) {
169 createBasicVector(entries, entryStride, weights, weightStrides);
200 int yStride = 1,
int zStride = 1,
bool usewgts =
false,
201 const scalar_t *wgts = NULL,
int wgtStride = 1)
202 : numIds_(numIds), idList_(ids), numEntriesPerID_(0),
203 numWeights_(usewgts == true) {
204 std::vector<const scalar_t *> values, weightValues;
205 std::vector<int> strides, weightStrides;
209 strides.push_back(xStride);
213 strides.push_back(yStride);
217 strides.push_back(zStride);
223 weightValues.push_back(wgts);
224 weightStrides.push_back(wgtStride);
226 createBasicVector(values, strides, weightValues, weightStrides);
242 auto hostIds = Kokkos::create_mirror_view(this->kIds_);
243 Kokkos::deep_copy(hostIds, this->kIds_);
260 "Invalid weight index.");
263 weights_[idx].getStridedList(length, weights, stride);
267 int idx = 0)
const override {
269 "Invalid weight index.");
272 typename Base::WeightsDeviceView1D(
"weights", kWeights_.extent(0));
275 hostWgts = Kokkos::create_mirror_view(weightsDevice);
276 Kokkos::deep_copy(hostWgts, weightsDevice);
280 auto hostWeights = Kokkos::create_mirror_view(kWeights_);
281 Kokkos::deep_copy(hostWeights, kWeights_);
286 int idx = 0)
const override {
288 "Invalid weight index.");
290 const auto size = kWeights_.extent(0);
291 deviceWgts =
typename Base::WeightsDeviceView1D(
"weights", size);
293 Kokkos::parallel_for(
294 size, KOKKOS_CLASS_LAMBDA(
const int id) {
295 deviceWgts(
id) = kWeights_(
id, idx);
314 if (idx < 0 || idx >= numEntriesPerID_) {
315 std::ostringstream emsg;
316 emsg << __FILE__ <<
":" << __LINE__ <<
" Invalid vector index " << idx
318 throw std::runtime_error(emsg.str());
321 entries_[idx].getStridedList(length, entries, stride);
330 &entries)
const override {
331 auto hostEntries = Kokkos::create_mirror_view(kEntries_);
332 Kokkos::deep_copy(hostEntries, kEntries_);
333 entries = hostEntries;
337 &entries)
const override {
343 const gno_t *idList_;
345 int numEntriesPerID_;
349 ArrayRCP<StridedData<lno_t, scalar_t>> entries_;
350 ArrayRCP<StridedData<lno_t, scalar_t>> weights_;
353 typename Base::IdsDeviceView kIds_;
355 typename Base::WeightsDeviceView kWeights_;
357 void createBasicVector(std::vector<const scalar_t *> &entries,
358 std::vector<int> &entryStride,
359 std::vector<const scalar_t *> &
weights,
360 std::vector<int> &weightStrides) {
365 kIds_ =
typename Base::IdsDeviceView(
366 Kokkos::ViewAllocateWithoutInitializing(
"ids"), numIds_);
367 auto host_kIds_ = Kokkos::create_mirror_view(kIds_);
368 for (
int n = 0; n < numIds_; ++n) {
369 host_kIds_(n) = idList_[n];
371 Kokkos::deep_copy(kIds_, host_kIds_);
375 entries_ = arcp(
new input_t[numEntriesPerID_], 0, numEntriesPerID_,
true);
376 for (
int v = 0; v < numEntriesPerID_; v++) {
377 if (entryStride.size())
378 stride = entryStride[v];
379 ArrayRCP<const scalar_t> eltV(entries[v], 0, stride * numIds_,
false);
380 entries_[v] = input_t(eltV, stride);
385 Kokkos::ViewAllocateWithoutInitializing(
"entries"), numIds_,
391 auto host_kokkos_entries = Kokkos::create_mirror_view(this->kEntries_);
393 for (
int idx = 0; idx < numEntriesPerID_; idx++) {
394 entries_[idx].getStridedList(length, entriesPtr, stride);
395 size_t fill_index = 0;
396 for (
int n = 0; n < numIds_; ++n) {
397 host_kokkos_entries(fill_index++, idx) = entriesPtr[n * stride];
400 Kokkos::deep_copy(this->kEntries_, host_kokkos_entries);
405 weights_ = arcp(
new input_t[numWeights_], 0, numWeights_,
true);
406 for (
int w = 0; w < numWeights_; w++) {
407 if (weightStrides.size())
408 stride = weightStrides[w];
409 ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride * numIds_,
false);
410 weights_[w] = input_t(wgtV, stride);
414 kWeights_ =
typename Base::WeightsDeviceView(
415 Kokkos::ViewAllocateWithoutInitializing(
"kokkos weights"), numIds_,
419 auto host_weight_temp_values =
420 Kokkos::create_mirror_view(this->kWeights_);
421 for (
int idx = 0; idx < numWeights_; ++idx) {
424 weights_[idx].getStridedList(length, weightsPtr, stride);
425 size_t fill_index = 0;
426 for (
size_t n = 0; n < length; n += stride) {
427 host_weight_temp_values(fill_index++, idx) = weightsPtr[n];
430 Kokkos::deep_copy(this->kWeights_, host_weight_temp_values);
void getEntriesDeviceView(typename AdapterWithCoords< User >::CoordsDeviceView &entries) const override
Provide a Kokkos view (Device side) to the elements of the specified vector.
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getIDsHostView(typename Base::ConstIdsHostView &ids) const override
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)
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.
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
map_t::global_ordinal_type gno_t
Defines the VectorAdapter interface.
void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
typename InputTraits< User >::part_t part_t
void getIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override
int getNumEntriesPerID() const
Return the number of vectors.
typename InputTraits< User >::node_t node_t
typename CoordsDeviceView::HostMirror CoordsHostView
void getWeightsHostView(typename Base::WeightsHostView &wgts) const override
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
typename InputTraits< User >::gno_t gno_t
void getEntriesKokkosView(typename AdapterWithCoords< User >::CoordsDeviceView &entries) const
VectorAdapter defines the interface for vector input.
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
The StridedData class manages lists of weights or coordinates.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
map_t::local_ordinal_type lno_t
void getEntriesHostView(typename AdapterWithCoords< User >::CoordsHostView &entries) const override
Provide a Kokkos view (Host side) to the elements of the specified vector.
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...
typename InputTraits< User >::offset_t offset_t
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
typename BaseAdapter< User >::scalar_t scalar_t
void getWeightsDeviceView(typename Base::WeightsDeviceView &wgts) const override
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
void getIDsKokkosView(typename Base::ConstIdsDeviceView &ids) const
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...
This file defines the StridedData class.