Zoltan2
Zoltan2_XpetraMultiVectorAdapter.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_XPETRAMULTIVECTORADAPTER_HPP_
51 #define _ZOLTAN2_XPETRAMULTIVECTORADAPTER_HPP_
52 
53 #include <Zoltan2_XpetraTraits.hpp>
55 #include <Zoltan2_StridedData.hpp>
57 
58 #include <Xpetra_EpetraMultiVector.hpp>
59 #include <Xpetra_TpetraMultiVector.hpp>
60 
61 namespace Zoltan2 {
62 
80 template <typename User>
81  class XpetraMultiVectorAdapter : public VectorAdapter<User> {
82 public:
83 
84 #ifndef DOXYGEN_SHOULD_SKIP_THIS
85  typedef typename InputTraits<User>::scalar_t scalar_t;
86  typedef typename InputTraits<User>::lno_t lno_t;
87  typedef typename InputTraits<User>::gno_t gno_t;
88  typedef typename InputTraits<User>::part_t part_t;
89  typedef typename InputTraits<User>::node_t node_t;
90  typedef VectorAdapter<User> base_adapter_t;
91  typedef User user_t;
92  typedef User userCoord_t;
93 
94  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
95  typedef Xpetra::TpetraMultiVector<
96  scalar_t, lno_t, gno_t, node_t> xt_mvector_t;
97  typedef Xpetra::EpetraMultiVector xe_mvector_t;
98 #endif
99 
103 
119  XpetraMultiVectorAdapter(const RCP<const User> &invector,
120  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides);
121 
127  XpetraMultiVectorAdapter(const RCP<const User> &invector);
128 
129 
131  // The Adapter interface.
133 
134  size_t getLocalNumIDs() const { return vector_->getLocalLength();}
135 
136  void getIDsView(const gno_t *&ids) const
137  {
138  ids = map_->getNodeElementList().getRawPtr();
139  }
140 
141  int getNumWeightsPerID() const { return numWeights_;}
142 
143  void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
144  {
145  env_->localInputAssertion(__FILE__, __LINE__, "invalid weight index",
146  idx >= 0 && idx < numWeights_, BASIC_ASSERTION);
147  size_t length;
148  weights_[idx].getStridedList(length, weights, stride);
149  }
150 
152  // The VectorAdapter interface.
154 
155  int getNumEntriesPerID() const {return vector_->getNumVectors();}
156 
157  void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const;
158 
159  template <typename Adapter>
160  void applyPartitioningSolution(const User &in, User *&out,
161  const PartitioningSolution<Adapter> &solution) const;
162 
163  template <typename Adapter>
164  void applyPartitioningSolution(const User &in, RCP<User> &out,
165  const PartitioningSolution<Adapter> &solution) const;
166 
167 private:
168 
169  RCP<const User> invector_;
170  RCP<const x_mvector_t> vector_;
171  RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
172  RCP<Environment> env_; // for error messages, etc.
173  lno_t base_;
174 
175  int numWeights_;
176  ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
177 };
178 
180 // Definitions
182 
183 template <typename User>
185  const RCP<const User> &invector,
186  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
187  invector_(invector), vector_(), map_(),
188  env_(rcp(new Environment)), base_(),
189  numWeights_(weights.size()), weights_(weights.size())
190 {
191  typedef StridedData<lno_t, scalar_t> input_t;
192 
193  RCP<x_mvector_t> tmp =
194  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
195  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
196  map_ = vector_->getMap();
197  base_ = map_->getIndexBase();
198 
199  size_t length = vector_->getLocalLength();
200 
201  if (length > 0 && numWeights_ > 0){
202  int stride = 1;
203  for (int w=0; w < numWeights_; w++){
204  if (weightStrides.size())
205  stride = weightStrides[w];
206  ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*length, false);
207  weights_[w] = input_t(wgtV, stride);
208  }
209  }
210 }
211 
212 
214 template <typename User>
216  const RCP<const User> &invector):
217  invector_(invector), vector_(), map_(),
218  env_(rcp(new Environment)), base_(),
219  numWeights_(0), weights_()
220 {
221  RCP<x_mvector_t> tmp =
222  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
223  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
224  map_ = vector_->getMap();
225  base_ = map_->getIndexBase();
226 }
227 
229 template <typename User>
231  const scalar_t *&elements, int &stride, int idx) const
232 {
233  size_t vecsize;
234  stride = 1;
235  elements = NULL;
236  if (map_->lib() == Xpetra::UseTpetra){
237  const xt_mvector_t *tvector =
238  dynamic_cast<const xt_mvector_t *>(vector_.get());
239 
240  vecsize = tvector->getLocalLength();
241  if (vecsize > 0){
242  ArrayRCP<const scalar_t> data = tvector->getData(idx);
243  elements = data.get();
244  }
245  }
246  else if (map_->lib() == Xpetra::UseEpetra){
247  const xe_mvector_t *evector =
248  dynamic_cast<const xe_mvector_t *>(vector_.get());
249 
250  vecsize = evector->getLocalLength();
251  if (vecsize > 0){
252  ArrayRCP<const double> data = evector->getData(idx);
253 
254  // Cast so this will compile when scalar_t is not double,
255  // a case when this code should never execute.
256  elements = reinterpret_cast<const scalar_t *>(data.get());
257  }
258  }
259  else{
260  throw std::logic_error("invalid underlying lib");
261  }
262 }
263 
265 template <typename User>
266  template <typename Adapter>
268  const User &in, User *&out,
269  const PartitioningSolution<Adapter> &solution) const
270 {
271  // Get an import list (rows to be received)
272  size_t numNewRows;
273  ArrayRCP<gno_t> importList;
274  try{
275  numNewRows = Zoltan2::getImportList<Adapter,
277  (solution, this, importList);
278  }
280 
281  // Move the rows, creating a new vector.
282  RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewRows,
283  importList.getRawPtr());
284  out = outPtr.get();
285  outPtr.release();
286 }
287 
289 template <typename User>
290  template <typename Adapter>
292  const User &in, RCP<User> &out,
293  const PartitioningSolution<Adapter> &solution) const
294 {
295  // Get an import list (rows to be received)
296  size_t numNewRows;
297  ArrayRCP<gno_t> importList;
298  try{
299  numNewRows = Zoltan2::getImportList<Adapter,
301  (solution, this, importList);
302  }
304 
305  // Move the rows, creating a new vector.
306  out = XpetraTraits<User>::doMigration(in, numNewRows,
307  importList.getRawPtr());
308 }
309 
310 } //namespace Zoltan2
311 
312 #endif
InputTraits< User >::scalar_t scalar_t
Helper functions for Partitioning Problems.
fast typical checks for valid arguments
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
InputTraits< User >::gno_t gno_t
default_part_t part_t
The data type to represent part numbers.
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor.
Defines the VectorAdapter interface.
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...
Traits of Xpetra classes, including migration method.
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
int getNumEntriesPerID() const
Return the number of vectors (typically one).
A PartitioningSolution is a solution to a partitioning problem.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) 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.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
InputTraits< User >::part_t part_t
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
An adapter for Xpetra::MultiVector.
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.
size_t getLocalNumIDs() const
Returns the number of objects on this process.
default_scalar_t scalar_t
The data type for weights and coordinates.
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
This file defines the StridedData class.