Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 
50 #ifndef _ZOLTAN2_BASICVECTORADAPTER_HPP_
51 #define _ZOLTAN2_BASICVECTORADAPTER_HPP_
52 
54 #include <Zoltan2_StridedData.hpp>
55 
56 namespace Zoltan2 {
57 
80 template <typename User>
81  class BasicVectorAdapter : public VectorAdapter<User> {
82 
83 public:
84 
85 #ifndef DOXYGEN_SHOULD_SKIP_THIS
86 
87  typedef typename InputTraits<User>::scalar_t scalar_t;
88  typedef typename InputTraits<User>::lno_t lno_t;
89  typedef typename InputTraits<User>::gno_t gno_t;
90  typedef typename InputTraits<User>::offset_t offset_t;
91  typedef typename InputTraits<User>::part_t part_t;
92  typedef typename InputTraits<User>::node_t node_t;
93  typedef User user_t;
94 
95 #endif
96 
112  BasicVectorAdapter(lno_t numIds, const gno_t *ids,
113  const scalar_t *entries, int entryStride=1,
114  bool usewgts=false,
115  const scalar_t *wgts=NULL, int wgtStride=1):
116  numIds_(numIds), idList_(ids),
117  numEntriesPerID_(1), entries_(),
118  numWeights_(usewgts==true), weights_()
119  {
120  std::vector<const scalar_t *> values;
121  std::vector<int> strides;
122  std::vector<const scalar_t *> weightValues;
123  std::vector<int> weightStrides;
124 
125  values.push_back(entries);
126  strides.push_back(entryStride);
127  if (usewgts) {
128  weightValues.push_back(wgts);
129  weightStrides.push_back(wgtStride);
130  }
131 
132  createBasicVector(values, strides, weightValues, weightStrides);
133  }
134 
160  BasicVectorAdapter(lno_t numIds, const gno_t *ids,
161  std::vector<const scalar_t *> &entries, std::vector<int> &entryStride,
162  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
163  numIds_(numIds), idList_(ids),
164  numEntriesPerID_(entries.size()), entries_(),
165  numWeights_(weights.size()), weights_()
166  {
167  createBasicVector(entries, entryStride, weights, weightStrides);
168  }
169 
196  BasicVectorAdapter(lno_t numIds, const gno_t *ids,
197  const scalar_t *x, const scalar_t *y,
198  const scalar_t *z,
199  int xStride=1, int yStride=1, int zStride=1,
200  bool usewgts=false, const scalar_t *wgts=NULL,
201  int wgtStride=1) :
202  numIds_(numIds), idList_(ids), numEntriesPerID_(0), entries_(),
203  numWeights_(usewgts==true), weights_()
204  {
205  std::vector<const scalar_t *> values, weightValues;
206  std::vector<int> strides, weightStrides;
207 
208  if (x){
209  values.push_back(x);
210  strides.push_back(xStride);
211  numEntriesPerID_++;
212  if (y){
213  values.push_back(y);
214  strides.push_back(yStride);
215  numEntriesPerID_++;
216  if (z){
217  values.push_back(z);
218  strides.push_back(zStride);
219  numEntriesPerID_++;
220  }
221  }
222  }
223  if (usewgts) {
224  weightValues.push_back(wgts);
225  weightStrides.push_back(wgtStride);
226  }
227  createBasicVector(values, strides, weightValues, weightStrides);
228  }
229 
230 
232 
234  // The Adapter interface.
236 
237  size_t getLocalNumIDs() const { return numIds_;}
238 
239  void getIDsView(const gno_t *&ids) const {ids = idList_;}
240 
241  int getNumWeightsPerID() const { return numWeights_;}
242 
243  void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
244  {
245  if (idx < 0 || idx >= numWeights_) {
246  std::ostringstream emsg;
247  emsg << __FILE__ << ":" << __LINE__
248  << " Invalid vector index " << idx << std::endl;
249  throw std::runtime_error(emsg.str());
250  }
251  size_t length;
252  weights_[idx].getStridedList(length, weights, stride);
253  }
254 
256  // The VectorAdapter interface.
258 
259  int getNumEntriesPerID() const { return numEntriesPerID_;}
260 
261  void getEntriesView(const scalar_t *&entries, int &stride, int idx = 0) const
262  {
263  if (idx < 0 || idx >= numEntriesPerID_) {
264  std::ostringstream emsg;
265  emsg << __FILE__ << ":" << __LINE__
266  << " Invalid vector index " << idx << std::endl;
267  throw std::runtime_error(emsg.str());
268  }
269  size_t length;
270  entries_[idx].getStridedList(length, entries, stride);
271  }
272 
273 private:
274 
275  lno_t numIds_;
276  const gno_t *idList_;
277 
278  int numEntriesPerID_;
279  ArrayRCP<StridedData<lno_t, scalar_t> > entries_ ;
280 
281  int numWeights_;
282  ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
283 
284  void createBasicVector(
285  std::vector<const scalar_t *> &entries, std::vector<int> &entryStride,
286  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides)
287  {
288  typedef StridedData<lno_t,scalar_t> input_t;
289 
290  if (numIds_){
291  int stride = 1;
292  entries_ = arcp(new input_t[numEntriesPerID_], 0, numEntriesPerID_, true);
293  for (int v=0; v < numEntriesPerID_; v++) {
294  if (entryStride.size()) stride = entryStride[v];
295  ArrayRCP<const scalar_t> eltV(entries[v], 0, stride*numIds_, false);
296  entries_[v] = input_t(eltV, stride);
297  }
298  }
299 
300  if (numWeights_) {
301  int stride = 1;
302  weights_ = arcp(new input_t [numWeights_], 0, numWeights_, true);
303  for (int w=0; w < numWeights_; w++){
304  if (weightStrides.size()) stride = weightStrides[w];
305  ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*numIds_, false);
306  weights_[w] = input_t(wgtV, stride);
307  }
308  }
309  }
310 };
311 
312 } //namespace Zoltan2
313 
314 #endif
InputTraits< User >::scalar_t scalar_t
size_t getLocalNumIDs() const
Returns the number of objects on this process.
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.
InputTraits< User >::gno_t gno_t
default_part_t part_t
The data type to represent part numbers.
default_offset_t offset_t
The data type to represent offsets.
InputTraits< User >::lno_t lno_t
Defines the VectorAdapter interface.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
int getNumEntriesPerID() const
Return the number of vectors.
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
default_lno_t lno_t
The ordinal type (e.g., int, long, int64_t) that represents local counts and local indices...
VectorAdapter defines the interface for vector input.
The StridedData class manages lists of weights or coordinates.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
InputTraits< User >::part_t part_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...
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.
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process&#39; identifiers.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
InputTraits< User >::offset_t offset_t
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...
default_scalar_t scalar_t
The data type for weights and coordinates.
This file defines the StridedData class.
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74