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 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
59 #include <Xpetra_EpetraMultiVector.hpp>
60 #endif
61 #include <Xpetra_TpetraMultiVector.hpp>
62 
63 namespace Zoltan2 {
64 
82 template <typename User>
83  class XpetraMultiVectorAdapter : public VectorAdapter<User> {
84 public:
85 
86 #ifndef DOXYGEN_SHOULD_SKIP_THIS
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>::part_t part_t;
91  typedef typename InputTraits<User>::node_t node_t;
92  typedef User user_t;
93  typedef User userCoord_t;
94 
95  typedef Xpetra::MultiVector<scalar_t, lno_t, gno_t, node_t> x_mvector_t;
96  typedef Xpetra::TpetraMultiVector<
97  scalar_t, lno_t, gno_t, node_t> xt_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 
142  Kokkos::View<const gno_t *, typename node_t::device_type> &ids) const {
143  if (map_->lib() == Xpetra::UseTpetra) {
144  using device_type = typename node_t::device_type;
145  const xt_mvector_t *tvector =
146  dynamic_cast<const xt_mvector_t *>(vector_.get());
147  // MJ can be running Host, CudaSpace, or CudaUVMSpace while Map now
148  // internally never stores CudaUVMSpace so we may need a conversion.
149  // However Map stores both Host and CudaSpace so this could be improved
150  // if device_type was CudaSpace. Then we could add a new accessor to
151  // Map such as getMyGlobalIndicesDevice() which could be direct assigned
152  // here. Since Tpetra is still UVM dependent that is not going to happen
153  // yet so just leaving this as Host to device_type conversion for now.
154  ids = Kokkos::create_mirror_view_and_copy(device_type(),
155  tvector->getTpetra_MultiVector()->getMap()->getMyGlobalIndices());
156  }
157  else if (map_->lib() == Xpetra::UseEpetra) {
158 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
159  // this will call getIDsView to get raw ptr and create a view from it
161 #else
162  throw std::logic_error("Epetra requested, but Trilinos is not "
163  "built with Epetra");
164 #endif
165  }
166  else {
167  throw std::logic_error("getIDsKokkosView called but not on Tpetra or Epetra!");
168  }
169  }
170 
171  int getNumWeightsPerID() const { return numWeights_;}
172 
173  void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
174  {
175  if(idx<0 || idx >= numWeights_)
176  {
177  std::ostringstream emsg;
178  emsg << __FILE__ << ":" << __LINE__
179  << " Invalid weight index " << idx << std::endl;
180  throw std::runtime_error(emsg.str());
181  }
182 
183  size_t length;
184  weights_[idx].getStridedList(length, weights, stride);
185  }
186 
187  void getWeightsKokkos2dView(Kokkos::View<scalar_t **,
188  typename node_t::device_type> &wgt) const {
189  typedef Kokkos::View<scalar_t**, typename node_t::device_type> view_t;
190  wgt = view_t("wgts", vector_->getLocalLength(), numWeights_);
191  typename view_t::HostMirror host_wgt = Kokkos::create_mirror_view(wgt);
192  for(int idx = 0; idx < numWeights_; ++idx) {
193  const scalar_t * weights;
194  size_t length;
195  int stride;
196  weights_[idx].getStridedList(length, weights, stride);
197  size_t fill_index = 0;
198  for(size_t n = 0; n < length; n += stride) {
199  host_wgt(fill_index++,idx) = weights[n];
200  }
201  }
202  Kokkos::deep_copy(wgt, host_wgt);
203  }
204 
206  // The VectorAdapter interface.
208 
209  int getNumEntriesPerID() const {return vector_->getNumVectors();}
210 
211  void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const;
212 
214  // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
215  Kokkos::View<scalar_t **, Kokkos::LayoutLeft,
216  typename node_t::device_type> & elements) const;
217 
218  template <typename Adapter>
219  void applyPartitioningSolution(const User &in, User *&out,
220  const PartitioningSolution<Adapter> &solution) const;
221 
222  template <typename Adapter>
223  void applyPartitioningSolution(const User &in, RCP<User> &out,
224  const PartitioningSolution<Adapter> &solution) const;
225 
226 private:
227 
228  RCP<const User> invector_;
229  RCP<const x_mvector_t> vector_;
230  RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
231 
232  int numWeights_;
233  ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
234 };
235 
237 // Definitions
239 
240 template <typename User>
242  const RCP<const User> &invector,
243  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
244  invector_(invector), vector_(), map_(),
245  numWeights_(weights.size()), weights_(weights.size())
246 {
247  typedef StridedData<lno_t, scalar_t> input_t;
248 
249  try {
250  RCP<x_mvector_t> tmp =
251  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
252  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
253  }
255 
256  map_ = vector_->getMap();
257 
258  size_t length = vector_->getLocalLength();
259 
260  if (length > 0 && numWeights_ > 0){
261  int stride = 1;
262  for (int w=0; w < numWeights_; w++){
263  if (weightStrides.size())
264  stride = weightStrides[w];
265  ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*length, false);
266  weights_[w] = input_t(wgtV, stride);
267  }
268  }
269 }
270 
271 
273 template <typename User>
275  const RCP<const User> &invector):
276  invector_(invector), vector_(), map_(),
277  numWeights_(0), weights_()
278 {
279  try {
280  RCP<x_mvector_t> tmp =
281  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
282  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
283  }
285 
286  map_ = vector_->getMap();
287 }
288 
290 template <typename User>
292  const scalar_t *&elements, int &stride, int idx) const
293 {
294  size_t vecsize;
295  stride = 1;
296  elements = NULL;
297  if (map_->lib() == Xpetra::UseTpetra){
298  const xt_mvector_t *tvector =
299  dynamic_cast<const xt_mvector_t *>(vector_.get());
300 
301  vecsize = tvector->getLocalLength();
302  if (vecsize > 0){
303  ArrayRCP<const scalar_t> data = tvector->getData(idx);
304  elements = data.get();
305  }
306  }
307  else if (map_->lib() == Xpetra::UseEpetra){
308 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
309  typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
310  const xe_mvector_t *evector =
311  dynamic_cast<const xe_mvector_t *>(vector_.get());
312 
313  vecsize = evector->getLocalLength();
314  if (vecsize > 0){
315  ArrayRCP<const double> data = evector->getData(idx);
316 
317  // Cast so this will compile when scalar_t is not double,
318  // a case when this code should never execute.
319  elements = reinterpret_cast<const scalar_t *>(data.get());
320  }
321 #else
322  throw std::logic_error("Epetra requested, but Trilinos is not "
323  "built with Epetra");
324 #endif
325  }
326  else{
327  throw std::logic_error("invalid underlying lib");
328  }
329 }
330 
332 template <typename User>
334  // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
335  Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> & elements) const
336 {
337  if (map_->lib() == Xpetra::UseTpetra){
338  const xt_mvector_t *tvector =
339  dynamic_cast<const xt_mvector_t *>(vector_.get());
340  // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
341  Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> view2d =
342  tvector->getTpetra_MultiVector()->template getLocalView<node_t>(Tpetra::Access::ReadWrite);
343  elements = view2d;
344  }
345  else if (map_->lib() == Xpetra::UseEpetra){
346 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
347  typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
348  const xe_mvector_t *evector =
349  dynamic_cast<const xe_mvector_t *>(vector_.get());
350  elements =
351  Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type>
352  ("elements", evector->getLocalLength(), evector->getNumVectors());
353  if(evector->getLocalLength() > 0) {
354  for(size_t idx = 0; idx < evector->getNumVectors(); ++idx) {
355  const scalar_t * ptr;
356  int stride;
357  getEntriesView(ptr, stride, idx);
358  for(size_t n = 0; n < evector->getLocalLength(); ++n) {
359  elements(n, idx) = ptr[n];
360  }
361  }
362  }
363 #else
364  throw std::logic_error("Epetra requested, but Trilinos is not "
365  "built with Epetra");
366 #endif
367  }
368  else {
369  throw std::logic_error("getEntriesKokkosView called but not using Tpetra or Epetra!");
370  }
371 }
372 
374 template <typename User>
375  template <typename Adapter>
377  const User &in, User *&out,
378  const PartitioningSolution<Adapter> &solution) const
379 {
380  // Get an import list (rows to be received)
381  size_t numNewRows;
382  ArrayRCP<gno_t> importList;
383  try{
384  numNewRows = Zoltan2::getImportList<Adapter,
386  (solution, this, importList);
387  }
389 
390  // Move the rows, creating a new vector.
391  RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewRows,
392  importList.getRawPtr());
393  out = outPtr.get();
394  outPtr.release();
395 }
396 
398 template <typename User>
399  template <typename Adapter>
401  const User &in, RCP<User> &out,
402  const PartitioningSolution<Adapter> &solution) const
403 {
404  // Get an import list (rows to be received)
405  size_t numNewRows;
406  ArrayRCP<gno_t> importList;
407  try{
408  numNewRows = Zoltan2::getImportList<Adapter,
410  (solution, this, importList);
411  }
413 
414  // Move the rows, creating a new vector.
415  out = XpetraTraits<User>::doMigration(in, numNewRows,
416  importList.getRawPtr());
417 }
418 
419 } //namespace Zoltan2
420 
421 #endif
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:74
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines the VectorAdapter interface.
Traits of Xpetra classes, including migration method.
InputTraits< User >::node_t node_t
InputTraits< User >::part_t part_t
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
virtual void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
Provide a Kokkos view to this process' identifiers.
InputTraits< User >::gno_t gno_t
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
void getWeightsKokkos2dView(Kokkos::View< scalar_t **, typename node_t::device_type > &wgt) const
void getEntriesKokkosView(Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &elements) const
Provide a Kokkos view to the elements of the specified vector.
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor.
void getIDsKokkosView(Kokkos::View< const gno_t *, typename node_t::device_type > &ids) const
int getNumEntriesPerID() const
Return the number of vectors.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Created by mbenlioglu on Aug 31, 2020.
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...
static ArrayRCP< ArrayRCP< zscalar_t > > weights
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.
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...
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.