Zoltan2
Loading...
Searching...
No Matches
Zoltan2_BasicVectorAdapter.hpp
Go to the documentation of this file.
1// @HEADER
2//
3// ***********************************************************************
4//
5// Zoltan2: A package of combinatorial algorithms for scientific computing
6// Copyright 2012 Sandia Corporation
7//
8// Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9// the U.S. Government retains certain rights in this software.
10//
11// Redistribution and use in source and binary forms, with or without
12// modification, are permitted provided that the following conditions are
13// met:
14//
15// 1. Redistributions of source code must retain the above copyright
16// notice, this list of conditions and the following disclaimer.
17//
18// 2. Redistributions in binary form must reproduce the above copyright
19// notice, this list of conditions and the following disclaimer in the
20// documentation and/or other materials provided with the distribution.
21//
22// 3. Neither the name of the Corporation nor the names of the
23// contributors may be used to endorse or promote products derived from
24// this software without specific prior written permission.
25//
26// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37//
38// Questions? Contact Karen Devine (kddevin@sandia.gov)
39// Erik Boman (egboman@sandia.gov)
40// Siva Rajamanickam (srajama@sandia.gov)
41//
42// ***********************************************************************
43//
44// @HEADER
45
49
50#ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
51#define _ZOLTAN2_BASICVECTORADAPTER_HPP_
52
55
56namespace Zoltan2 {
57
84
85template <typename User>
86 class BasicVectorAdapter : public VectorAdapter<User> {
87
88public:
89
90#ifndef DOXYGEN_SHOULD_SKIP_THIS
91
93 typedef typename InputTraits<User>::lno_t lno_t;
94 typedef typename InputTraits<User>::gno_t gno_t;
96 typedef typename InputTraits<User>::part_t part_t;
97 typedef typename InputTraits<User>::node_t node_t;
98 typedef User user_t;
99
100#endif
101
116
117 BasicVectorAdapter(lno_t numIds, const gno_t *ids,
118 const scalar_t *entries, int entryStride=1,
119 bool usewgts=false,
120 const scalar_t *wgts=NULL, int wgtStride=1):
121 numIds_(numIds), idList_(ids),
122 numEntriesPerID_(1), entries_(),
123 numWeights_(usewgts==true), weights_()
124 {
125 std::vector<const scalar_t *> values;
126 std::vector<int> strides;
127 std::vector<const scalar_t *> weightValues;
128 std::vector<int> weightStrides;
129
130 values.push_back(entries);
131 strides.push_back(entryStride);
132 if (usewgts) {
133 weightValues.push_back(wgts);
134 weightStrides.push_back(wgtStride);
135 }
136
137 createBasicVector(values, strides, weightValues, weightStrides);
138 }
139
164
165 BasicVectorAdapter(lno_t numIds, const gno_t *ids,
166 std::vector<const scalar_t *> &entries, std::vector<int> &entryStride,
167 std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
168 numIds_(numIds), idList_(ids),
169 numEntriesPerID_(entries.size()), entries_(),
170 numWeights_(weights.size()), weights_()
171 {
172 createBasicVector(entries, entryStride, weights, weightStrides);
173 }
174
200
201 BasicVectorAdapter(lno_t numIds, const gno_t *ids,
202 const scalar_t *x, const scalar_t *y,
203 const scalar_t *z,
204 int xStride=1, int yStride=1, int zStride=1,
205 bool usewgts=false, const scalar_t *wgts=NULL,
206 int wgtStride=1) :
207 numIds_(numIds), idList_(ids), numEntriesPerID_(0), entries_(),
208 numWeights_(usewgts==true), weights_()
209 {
210 std::vector<const scalar_t *> values, weightValues;
211 std::vector<int> strides, weightStrides;
212
213 if (x){
214 values.push_back(x);
215 strides.push_back(xStride);
216 numEntriesPerID_++;
217 if (y){
218 values.push_back(y);
219 strides.push_back(yStride);
220 numEntriesPerID_++;
221 if (z){
222 values.push_back(z);
223 strides.push_back(zStride);
224 numEntriesPerID_++;
225 }
226 }
227 }
228 if (usewgts) {
229 weightValues.push_back(wgts);
230 weightStrides.push_back(wgtStride);
231 }
232 createBasicVector(values, strides, weightValues, weightStrides);
233 }
234
235
237
239 // The Adapter interface.
241
242 size_t getLocalNumIDs() const { return numIds_;}
243
244 void getIDsView(const gno_t *&ids) const {ids = idList_;}
245
246 void getIDsKokkosView(Kokkos::View<const gno_t *,
247 typename node_t::device_type> &ids) const
248 {
249 ids = this->kokkos_ids_;
250 }
251
252 int getNumWeightsPerID() const { return numWeights_;}
253
254 virtual void getWeightsKokkos2dView(Kokkos::View<scalar_t **,
255 typename node_t::device_type> &wgt) const
256 {
257 wgt = kokkos_weights_;
258 }
259
260 void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
261 {
262 if (idx < 0 || idx >= numWeights_) {
263 std::ostringstream emsg;
264 emsg << __FILE__ << ":" << __LINE__
265 << " Invalid vector index " << idx << std::endl;
266 throw std::runtime_error(emsg.str());
267 }
268 size_t length;
269 weights_[idx].getStridedList(length, weights, stride);
270 }
271
273 // The VectorAdapter interface.
275
276 int getNumEntriesPerID() const { return numEntriesPerID_;}
277
278 void getEntriesView(const scalar_t *&entries, int &stride, int idx = 0) const
279 {
280 if (idx < 0 || idx >= numEntriesPerID_) {
281 std::ostringstream emsg;
282 emsg << __FILE__ << ":" << __LINE__
283 << " Invalid vector index " << idx << std::endl;
284 throw std::runtime_error(emsg.str());
285 }
286 size_t length;
287 entries_[idx].getStridedList(length, entries, stride);
288 }
289
291 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
292 Kokkos::View<scalar_t **, Kokkos::LayoutLeft,
293 typename node_t::device_type> & entries) const
294 {
295 entries = kokkos_entries_;
296 }
297
298private:
299 lno_t numIds_;
300 const gno_t *idList_;
301
302 int numEntriesPerID_;
303 ArrayRCP<StridedData<lno_t, scalar_t> > entries_ ;
304
305 Kokkos::View<gno_t *, typename node_t::device_type> kokkos_ids_;
306
307 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
308 Kokkos::View<scalar_t **, Kokkos::LayoutLeft,
309 typename node_t::device_type> kokkos_entries_;
310
311 int numWeights_;
312 ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
313
314 Kokkos::View<scalar_t**, typename node_t::device_type> kokkos_weights_;
315
316 void createBasicVector(
317 std::vector<const scalar_t *> &entries, std::vector<int> &entryStride,
318 std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides)
319 {
320 typedef StridedData<lno_t,scalar_t> input_t;
321
322 if (numIds_){
323 // make kokkos ids
324 kokkos_ids_ = Kokkos::View<gno_t *, typename node_t::device_type>(
325 Kokkos::ViewAllocateWithoutInitializing("ids"), numIds_);
326 auto host_kokkos_ids_ = Kokkos::create_mirror_view(kokkos_ids_);
327 for(int n = 0; n < numIds_; ++n) {
328 host_kokkos_ids_(n) = idList_[n];
329 }
330 Kokkos::deep_copy(kokkos_ids_, host_kokkos_ids_);
331
332 // make coordinates
333 int stride = 1;
334 entries_ = arcp(new input_t[numEntriesPerID_], 0, numEntriesPerID_, true);
335 for (int v=0; v < numEntriesPerID_; v++) {
336 if (entryStride.size()) stride = entryStride[v];
337 ArrayRCP<const scalar_t> eltV(entries[v], 0, stride*numIds_, false);
338 entries_[v] = input_t(eltV, stride);
339 }
340
341 // setup kokkos entries
342 // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
343 kokkos_entries_ = Kokkos::View<scalar_t **, Kokkos::LayoutLeft,
344 typename node_t::device_type>(
345 Kokkos::ViewAllocateWithoutInitializing("entries"),
346 numIds_, numEntriesPerID_);
347
348 size_t length;
349 const scalar_t * entriesPtr;
350
351 auto host_kokkos_entries =
352 Kokkos::create_mirror_view(this->kokkos_entries_);
353
354 for (int idx=0; idx < numEntriesPerID_; idx++) {
355 entries_[idx].getStridedList(length, entriesPtr, stride);
356 size_t fill_index = 0;
357 for(int n = 0; n < numIds_; ++n) {
358 host_kokkos_entries(fill_index++,idx) = entriesPtr[n*stride];
359 }
360 }
361 Kokkos::deep_copy(this->kokkos_entries_, host_kokkos_entries);
362 }
363
364 if (numWeights_) {
365 int stride = 1;
366 weights_ = arcp(new input_t [numWeights_], 0, numWeights_, true);
367 for (int w=0; w < numWeights_; w++){
368 if (weightStrides.size()) stride = weightStrides[w];
369 ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*numIds_, false);
370 weights_[w] = input_t(wgtV, stride);
371 }
372
373 // set up final view with weights
374 kokkos_weights_ = Kokkos::View<scalar_t**,
375 typename node_t::device_type>(
376 Kokkos::ViewAllocateWithoutInitializing("kokkos weights"),
377 numIds_, numWeights_);
378
379 // setup weights
380 auto host_weight_temp_values =
381 Kokkos::create_mirror_view(this->kokkos_weights_);
382 for(int idx = 0; idx < numWeights_; ++idx) {
383 const scalar_t * weightsPtr;
384 size_t length;
385 weights_[idx].getStridedList(length, weightsPtr, stride);
386 size_t fill_index = 0;
387 for(size_t n = 0; n < length; n += stride) {
388 host_weight_temp_values(fill_index++,idx) = weightsPtr[n];
389 }
390 }
391 Kokkos::deep_copy(this->kokkos_weights_, host_weight_temp_values);
392 }
393 }
394};
395
396} //namespace Zoltan2
397
398#endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition Metric.cpp:74
This file defines the StridedData class.
Defines the VectorAdapter interface.
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...
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
virtual void getWeightsKokkos2dView(Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getEntriesKokkosView(Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &entries) const
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
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 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.