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 
141  int getNumWeightsPerID() const { return numWeights_;}
142 
143  void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
144  {
145  if(idx<0 || idx >= numWeights_)
146  {
147  std::ostringstream emsg;
148  emsg << __FILE__ << ":" << __LINE__
149  << " Invalid weight index " << idx << std::endl;
150  throw std::runtime_error(emsg.str());
151  }
152 
153  size_t length;
154  weights_[idx].getStridedList(length, weights, stride);
155  }
156 
158  // The VectorAdapter interface.
160 
161  int getNumEntriesPerID() const {return vector_->getNumVectors();}
162 
163  void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const;
164 
165  template <typename Adapter>
166  void applyPartitioningSolution(const User &in, User *&out,
167  const PartitioningSolution<Adapter> &solution) const;
168 
169  template <typename Adapter>
170  void applyPartitioningSolution(const User &in, RCP<User> &out,
171  const PartitioningSolution<Adapter> &solution) const;
172 
173 private:
174 
175  RCP<const User> invector_;
176  RCP<const x_mvector_t> vector_;
177  RCP<const Xpetra::Map<lno_t, gno_t, node_t> > map_;
178 
179  int numWeights_;
180  ArrayRCP<StridedData<lno_t, scalar_t> > weights_;
181 };
182 
184 // Definitions
186 
187 template <typename User>
189  const RCP<const User> &invector,
190  std::vector<const scalar_t *> &weights, std::vector<int> &weightStrides):
191  invector_(invector), vector_(), map_(),
192  numWeights_(weights.size()), weights_(weights.size())
193 {
194  typedef StridedData<lno_t, scalar_t> input_t;
195 
196  try {
197  RCP<x_mvector_t> tmp =
198  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
199  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
200  }
202 
203  map_ = vector_->getMap();
204 
205  size_t length = vector_->getLocalLength();
206 
207  if (length > 0 && numWeights_ > 0){
208  int stride = 1;
209  for (int w=0; w < numWeights_; w++){
210  if (weightStrides.size())
211  stride = weightStrides[w];
212  ArrayRCP<const scalar_t> wgtV(weights[w], 0, stride*length, false);
213  weights_[w] = input_t(wgtV, stride);
214  }
215  }
216 }
217 
218 
220 template <typename User>
222  const RCP<const User> &invector):
223  invector_(invector), vector_(), map_(),
224  numWeights_(0), weights_()
225 {
226  try {
227  RCP<x_mvector_t> tmp =
228  XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(invector));
229  vector_ = rcp_const_cast<const x_mvector_t>(tmp);
230  }
232 
233  map_ = vector_->getMap();
234 }
235 
237 template <typename User>
239  const scalar_t *&elements, int &stride, int idx) const
240 {
241  size_t vecsize;
242  stride = 1;
243  elements = NULL;
244  if (map_->lib() == Xpetra::UseTpetra){
245  const xt_mvector_t *tvector =
246  dynamic_cast<const xt_mvector_t *>(vector_.get());
247 
248  vecsize = tvector->getLocalLength();
249  if (vecsize > 0){
250  ArrayRCP<const scalar_t> data = tvector->getData(idx);
251  elements = data.get();
252  }
253  }
254  else if (map_->lib() == Xpetra::UseEpetra){
255 #if defined(HAVE_ZOLTAN2_EPETRA) && defined(HAVE_XPETRA_EPETRA)
256  typedef Xpetra::EpetraMultiVectorT<gno_t,node_t> xe_mvector_t;
257  const xe_mvector_t *evector =
258  dynamic_cast<const xe_mvector_t *>(vector_.get());
259 
260  vecsize = evector->getLocalLength();
261  if (vecsize > 0){
262  ArrayRCP<const double> data = evector->getData(idx);
263 
264  // Cast so this will compile when scalar_t is not double,
265  // a case when this code should never execute.
266  elements = reinterpret_cast<const scalar_t *>(data.get());
267  }
268 #else
269  throw std::logic_error("Epetra requested, but Trilinos is not "
270  "built with Epetra");
271 #endif
272  }
273  else{
274  throw std::logic_error("invalid underlying lib");
275  }
276 }
277 
279 template <typename User>
280  template <typename Adapter>
282  const User &in, User *&out,
283  const PartitioningSolution<Adapter> &solution) const
284 {
285  // Get an import list (rows to be received)
286  size_t numNewRows;
287  ArrayRCP<gno_t> importList;
288  try{
289  numNewRows = Zoltan2::getImportList<Adapter,
291  (solution, this, importList);
292  }
294 
295  // Move the rows, creating a new vector.
296  RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewRows,
297  importList.getRawPtr());
298  out = outPtr.get();
299  outPtr.release();
300 }
301 
303 template <typename User>
304  template <typename Adapter>
306  const User &in, RCP<User> &out,
307  const PartitioningSolution<Adapter> &solution) const
308 {
309  // Get an import list (rows to be received)
310  size_t numNewRows;
311  ArrayRCP<gno_t> importList;
312  try{
313  numNewRows = Zoltan2::getImportList<Adapter,
315  (solution, this, importList);
316  }
318 
319  // Move the rows, creating a new vector.
320  out = XpetraTraits<User>::doMigration(in, numNewRows,
321  importList.getRawPtr());
322 }
323 
324 } //namespace Zoltan2
325 
326 #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 >::part_t part_t
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
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 getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
XpetraMultiVectorAdapter(const RCP< const User > &invector, std::vector< const scalar_t * > &weights, std::vector< int > &weightStrides)
Constructor.
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.
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.