Zoltan2
Zoltan2_CoordinateModel.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 
51 #ifndef _ZOLTAN2_COORDINATEMODEL_HPP_
52 #define _ZOLTAN2_COORDINATEMODEL_HPP_
53 
54 #include <Zoltan2_Model.hpp>
55 #include <Zoltan2_MeshAdapter.hpp>
57 #include <Zoltan2_GraphAdapter.hpp>
60 #include <Zoltan2_StridedData.hpp>
61 
62 namespace Zoltan2 {
63 
70 template <typename Adapter>
71 class CoordinateModel : public Model<Adapter>
72 {
73 public:
74 
75 #ifndef DOXYGEN_SHOULD_SKIP_THIS
76  typedef typename Adapter::scalar_t scalar_t;
77  typedef typename Adapter::gno_t gno_t;
78  typedef typename Adapter::lno_t lno_t;
79  typedef typename Adapter::node_t node_t;
80  typedef typename Adapter::user_t user_t;
81  typedef typename Adapter::userCoord_t userCoord_t;
82  typedef StridedData<lno_t, scalar_t> input_t;
83 #endif
84 
86  // Constructors for each Adapter type
88 
89  // VectorAdapter
90  CoordinateModel(const RCP<const VectorAdapter<user_t> > &ia,
91  const RCP<const Environment> &env,
92  const RCP<const Comm<int> > &comm,
93  modelFlag_t &flags):
94  numGlobalCoordinates_(), env_(env), comm_(comm),
95  coordinateDim_(), gids_(),
96  xyz_(), userNumWeights_(0), weights_()
97  {
98  typedef VectorAdapter<user_t> adapterWithCoords_t;
99  sharedConstructor<adapterWithCoords_t>(&(*ia), env, comm, flags);
100  }
101 
102  // MatrixAdapter
104  const RCP<const Environment> &env,
105  const RCP<const Comm<int> > &comm,
106  modelFlag_t &flags) :
107  numGlobalCoordinates_(), env_(env), comm_(comm),
108  coordinateDim_(), gids_(),
109  xyz_(), userNumWeights_(0), weights_()
110  {
111  if (!(ia->coordinatesAvailable()))
112  throw std::logic_error("No coordinate info provided to MatrixAdapter.");
113  else {
114  typedef VectorAdapter<userCoord_t> adapterWithCoords_t;
115  adapterWithCoords_t *va = ia->getCoordinateInput();
116  sharedConstructor<adapterWithCoords_t>(va, env, comm, flags);
117  }
118  }
119 
120  // GraphAdapter
122  const RCP<const Environment> &env,
123  const RCP<const Comm<int> > &comm,
124  modelFlag_t &flags) :
125  numGlobalCoordinates_(), env_(env), comm_(comm),
126  coordinateDim_(), gids_(),
127  xyz_(), userNumWeights_(0), weights_()
128  {
129  if (!(ia->coordinatesAvailable()))
130  throw std::logic_error("No coordinate info provided to GraphAdapter.");
131  else {
132  typedef VectorAdapter<userCoord_t> adapterWithCoords_t;
133  adapterWithCoords_t *va = ia->getCoordinateInput();
134  sharedConstructor<adapterWithCoords_t>(va, env, comm, flags);
135  }
136  }
137 
138  // MeshAdapter
139  CoordinateModel(const RCP<const MeshAdapter<user_t> > &ia,
140  const RCP<const Environment> &env,
141  const RCP<const Comm<int> > &comm,
142  modelFlag_t &flags) :
143  numGlobalCoordinates_(), env_(env), comm_(comm),
144  coordinateDim_(), gids_(),
145  xyz_(), userNumWeights_(0), weights_()
146  {
147  typedef MeshAdapter<user_t> adapterWithCoords_t;
148  sharedConstructor<adapterWithCoords_t>(&(*ia), env, comm, flags);
149  }
150 
151  // IdentifierAdapter
153  const RCP<const Environment> &env,
154  const RCP<const Comm<int> > &comm,
155  modelFlag_t &flags)
156  {
157  throw std::logic_error(
158  "A coordinate model can not be build from an IdentifierAdapter");
159  }
160 
162  // CoordinateModel interface.
164 
167  int getCoordinateDim() const { return coordinateDim_;}
168 
171  size_t getLocalNumCoordinates() const { return gids_.size();}
172 
175  global_size_t getGlobalNumCoordinates() const {return numGlobalCoordinates_;}
176 
179  int getNumWeightsPerCoordinate() const { return userNumWeights_;}
180 
203  size_t getCoordinates(ArrayView<const gno_t> &Ids,
204  ArrayView<input_t> &xyz,
205  ArrayView<input_t> &wgts) const
206  {
207  xyz = xyz_.view(0, coordinateDim_);
208  wgts = weights_.view(0, userNumWeights_);
209 
210  size_t nCoord = getLocalNumCoordinates();
211  Ids = ArrayView<const gno_t>();
212 
213  if (nCoord){
214  Ids = Teuchos::arrayView<const gno_t>(
215  reinterpret_cast<const gno_t *>(gids_.getRawPtr()), nCoord);
216  }
217 
218  return nCoord;
219  }
220 
229  Kokkos::View<const gno_t *, typename node_t::device_type> &Ids,
230  // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
231  Kokkos::View<scalar_t **,
232  Kokkos::LayoutLeft, typename node_t::device_type> &xyz,
233  Kokkos::View<scalar_t **, typename node_t::device_type> &wgts) const
234  {
235  Ids = kokkos_gids_;
236  xyz = kokkos_xyz_;
237  wgts = kokkos_weights_;
238  return getLocalNumCoordinates();
239  }
240 
242  // The Model interface.
244 
245  size_t getLocalNumObjects() const
246  {
247  return getLocalNumCoordinates();
248  }
249 
250  size_t getGlobalNumObjects() const
251  {
252  return getGlobalNumCoordinates();
253  }
254 
255 private:
256  size_t numGlobalCoordinates_;
257  const RCP<const Environment> env_;
258  const RCP<const Comm<int> > comm_;
259  int coordinateDim_;
260 
261  // TODO: We now have a Kokkos version and non kokkos version so need to clean
262  // this up and perhaps eliminate the non-kokkos version completely.
263  // However not all tests are converted to Kokkos so keeping both forms around
264  // for now is probably necessary.
265  Kokkos::View<gno_t *, typename node_t::device_type> kokkos_gids_;
266  // coordinates in MJ are LayoutLeft since Tpetra Multivector gives LayoutLeft
267  Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> kokkos_xyz_;
268  Kokkos::View<scalar_t **, typename node_t::device_type> kokkos_weights_;
269 
270  ArrayRCP<const gno_t> gids_;
271  ArrayRCP<input_t> xyz_;
272  int userNumWeights_;
273  ArrayRCP<input_t> weights_;
274 
275  template <typename AdapterWithCoords>
276  void sharedConstructor(const AdapterWithCoords *ia,
277  const RCP<const Environment> &env,
278  const RCP<const Comm<int> > &comm,
279  modelFlag_t &flags);
280 
281 };
282 
283 
285 
286 // sharedConstructor
287 template <typename Adapter>
288 template <typename AdapterWithCoords>
289 void CoordinateModel<Adapter>::sharedConstructor(
290  const AdapterWithCoords *ia,
291  const RCP<const Environment> &/* env */,
292  const RCP<const Comm<int> > &comm,
293  modelFlag_t &/* flags */)
294 {
295  size_t nLocalIds = ia->getLocalNumIDs();
296 
297  // Get coordinates and weights (if any)
298 
299  int tmp[2], gtmp[2];
300  tmp[0] = ia->getDimension();
301  tmp[1] = ia->getNumWeightsPerID();
302  Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 2, tmp, gtmp);
303  coordinateDim_ = gtmp[0];
304  userNumWeights_ = gtmp[1];
305 
306  env_->localBugAssertion(__FILE__, __LINE__, "coordinate dimension",
307  coordinateDim_ > 0, COMPLEX_ASSERTION);
308 
309  input_t *coordArray = new input_t [coordinateDim_];
310  input_t *weightArray = NULL;
311  if (userNumWeights_)
312  weightArray = new input_t [userNumWeights_];
313 
314  env_->localMemoryAssertion(__FILE__, __LINE__, userNumWeights_+coordinateDim_,
315  coordArray && (!userNumWeights_|| weightArray));
316 
317 
318  if (nLocalIds){
319  const gno_t *gids=NULL;
320 
321  ia->getIDsView(gids);
322  gids_ = arcp(gids, 0, nLocalIds, false);
323 
324  for (int dim=0; dim < coordinateDim_; dim++){
325  int stride;
326  const scalar_t *coords=NULL;
327  try{
328  ia->getCoordinatesView(coords, stride, dim);
329  }
331 
332  ArrayRCP<const scalar_t> cArray(coords, 0, nLocalIds*stride, false);
333  coordArray[dim] = input_t(cArray, stride);
334  }
335 
336  for (int idx=0; idx < userNumWeights_; idx++){
337  int stride;
338  const scalar_t *weights;
339  try{
340  ia->getWeightsView(weights, stride, idx);
341  }
343 
344  ArrayRCP<const scalar_t> wArray(weights, 0, nLocalIds*stride, false);
345  weightArray[idx] = input_t(wArray, stride);
346  }
347  }
348 
349  xyz_ = arcp(coordArray, 0, coordinateDim_);
350 
351  if (userNumWeights_)
352  weights_ = arcp(weightArray, 0, userNumWeights_);
353 
354  // These are deep copies so we don't hold on to refs of the device views causing problems without UVM
355  Kokkos::View<const gno_t *, typename node_t::device_type> kokkos_gids;
356  ia->getIDsKokkosView(kokkos_gids);
357  kokkos_gids_ = Kokkos::View<gno_t *, typename node_t::device_type>("kokkos_gids_",kokkos_gids.extent(0));
358  Kokkos::deep_copy(kokkos_gids_, kokkos_gids);
359 
360  Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type> kokkos_xyz;
361  ia->getCoordinatesKokkosView(kokkos_xyz);
362  kokkos_xyz_ = Kokkos::View<scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type>("kokkos_xyz", kokkos_xyz.extent(0), kokkos_xyz.extent(1));
363  Kokkos::deep_copy(kokkos_xyz_, kokkos_xyz);
364 
365  if(userNumWeights_ > 0) {
366  Kokkos::View<scalar_t **, typename node_t::device_type> kokkos_weights;
367  ia->getWeightsKokkosView(kokkos_weights);
368  kokkos_weights_ = Kokkos::View<scalar_t **, typename node_t::device_type>("kokkos_weights_",kokkos_weights.extent(0), kokkos_weights.extent(1));
369  Kokkos::deep_copy(kokkos_weights_, kokkos_weights);
370  }
371 
372 
373  Teuchos::reduceAll<int, size_t>(*comm, Teuchos::REDUCE_SUM, 1,
374  &nLocalIds, &numGlobalCoordinates_);
375 
376  env_->memory("After construction of coordinate model");
377 }
378 
379 } // namespace Zoltan2
380 
381 #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.
Defines the GraphAdapter interface.
Defines the IdentifierAdapter interface.
Defines the MatrixAdapter interface.
Defines the MeshAdapter interface.
Defines the Model interface.
This file defines the StridedData class.
Defines the VectorAdapter interface.
This class provides geometric coordinates with optional weights to the Zoltan2 algorithm.
global_size_t getGlobalNumCoordinates() const
Returns the global number coordinates.
size_t getGlobalNumObjects() const
Return the global number of objects.
CoordinateModel(const RCP< const VectorAdapter< user_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
size_t getLocalNumObjects() const
Return the local number of objects.
CoordinateModel(const RCP< const IdentifierAdapter< user_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
CoordinateModel(const RCP< const MeshAdapter< user_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
size_t getCoordinatesKokkos(Kokkos::View< const gno_t *, typename node_t::device_type > &Ids, Kokkos::View< scalar_t **, Kokkos::LayoutLeft, typename node_t::device_type > &xyz, Kokkos::View< scalar_t **, typename node_t::device_type > &wgts) const
Returns the coordinate ids, values and optional weights.
int getCoordinateDim() const
Returns the dimension of the coordinates.
CoordinateModel(const RCP< const GraphAdapter< user_t, userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
CoordinateModel(const RCP< const MatrixAdapter< user_t, userCoord_t > > &ia, const RCP< const Environment > &env, const RCP< const Comm< int > > &comm, modelFlag_t &flags)
size_t getCoordinates(ArrayView< const gno_t > &Ids, ArrayView< input_t > &xyz, ArrayView< input_t > &wgts) const
Returns the coordinate ids, values and optional weights.
size_t getLocalNumCoordinates() const
Returns the number of coordinates on this process.
int getNumWeightsPerCoordinate() const
Returns the number (0 or greater) of weights per coordinate.
GraphAdapter defines the interface for graph-based user data.
IdentifierAdapter defines the interface for identifiers.
MatrixAdapter defines the adapter interface for matrices.
MeshAdapter defines the interface for mesh input.
The base class for all model classes.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Created by mbenlioglu on Aug 31, 2020.
Tpetra::global_size_t global_size_t
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
@ COMPLEX_ASSERTION
more involved, like validate a graph
static ArrayRCP< ArrayRCP< zscalar_t > > weights