Zoltan2
Zoltan2_TpetraRowMatrixAdapter.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_TPETRAROWMATRIXADAPTER_HPP_
51 #define _ZOLTAN2_TPETRAROWMATRIXADAPTER_HPP_
52 
54 #include <Zoltan2_StridedData.hpp>
56 
57 #include <Tpetra_RowMatrix.hpp>
58 
59 namespace Zoltan2 {
60 
62 
75 template <typename User, typename UserCoord=User>
76  class TpetraRowMatrixAdapter : public MatrixAdapter<User,UserCoord> {
77 public:
78 
79 #ifndef DOXYGEN_SHOULD_SKIP_THIS
80  typedef typename InputTraits<User>::scalar_t scalar_t;
81  typedef typename InputTraits<User>::offset_t offset_t;
82  typedef typename InputTraits<User>::lno_t lno_t;
83  typedef typename InputTraits<User>::gno_t gno_t;
84  typedef typename InputTraits<User>::part_t part_t;
85  typedef typename InputTraits<User>::node_t node_t;
86  typedef User user_t;
87  typedef UserCoord userCoord_t;
88 #endif
89 
93 
99  TpetraRowMatrixAdapter(const RCP<const User> &inmatrix,
100  int nWeightsPerRow=0);
101 
114  void setWeights(const scalar_t *weightVal, int stride, int idx = 0);
115 
131  void setRowWeights(const scalar_t *weightVal, int stride, int idx = 0);
132 
138  void setWeightIsDegree(int idx);
139 
145  void setRowWeightIsNumberOfNonZeros(int idx);
146 
148  // The MatrixAdapter interface.
150 
151  size_t getLocalNumRows() const {
152  return matrix_->getNodeNumRows();
153  }
154 
155  size_t getLocalNumColumns() const {
156  return matrix_->getNodeNumCols();
157  }
158 
159  size_t getLocalNumEntries() const {
160  return matrix_->getNodeNumEntries();
161  }
162 
163  bool CRSViewAvailable() const { return true; }
164 
165  void getRowIDsView(const gno_t *&rowIds) const
166  {
167  ArrayView<const gno_t> rowView = rowMap_->getNodeElementList();
168  rowIds = rowView.getRawPtr();
169  }
170 
171  void getCRSView(ArrayRCP<const offset_t> &offsets, ArrayRCP<const gno_t> &colIds) const
172  {
173  offsets = offset_;
174  colIds = columnIds_;
175  }
176 
177  void getCRSView(ArrayRCP<const offset_t> &offsets,
178  ArrayRCP<const gno_t> &colIds,
179  ArrayRCP<const scalar_t> &values) const
180  {
181  offsets = offset_;
182  colIds = columnIds_;
183  values = values_;
184  }
185 
186 
187  int getNumWeightsPerRow() const { return nWeightsPerRow_; }
188 
189  void getRowWeightsView(const scalar_t *&weights, int &stride,
190  int idx = 0) const
191  {
192  if(idx<0 || idx >= nWeightsPerRow_)
193  {
194  std::ostringstream emsg;
195  emsg << __FILE__ << ":" << __LINE__
196  << " Invalid row weight index " << idx << std::endl;
197  throw std::runtime_error(emsg.str());
198  }
199 
200 
201  size_t length;
202  rowWeights_[idx].getStridedList(length, weights, stride);
203  }
204 
205  bool useNumNonzerosAsRowWeight(int idx) const { return numNzWeight_[idx];}
206 
207  template <typename Adapter>
208  void applyPartitioningSolution(const User &in, User *&out,
209  const PartitioningSolution<Adapter> &solution) const;
210 
211  template <typename Adapter>
212  void applyPartitioningSolution(const User &in, RCP<User> &out,
213  const PartitioningSolution<Adapter> &solution) const;
214 
215 private:
216 
217  RCP<const User> matrix_;
218  RCP<const Tpetra::Map<lno_t, gno_t, node_t> > rowMap_;
219  RCP<const Tpetra::Map<lno_t, gno_t, node_t> > colMap_;
220  ArrayRCP<offset_t> offset_;
221  ArrayRCP<gno_t> columnIds_; // TODO: KDD Is it necessary to copy and store
222  ArrayRCP<scalar_t> values_; // TODO: the matrix here? Would prefer views.
223 
224  int nWeightsPerRow_;
225  ArrayRCP<StridedData<lno_t, scalar_t> > rowWeights_;
226  ArrayRCP<bool> numNzWeight_;
227 
228  bool mayHaveDiagonalEntries;
229 
230  RCP<User> doMigration(const User &from, size_t numLocalRows,
231  const gno_t *myNewRows) const;
232 };
233 
235 // Definitions
237 
238 template <typename User, typename UserCoord>
240  const RCP<const User> &inmatrix, int nWeightsPerRow):
241  matrix_(inmatrix), rowMap_(), colMap_(),
242  offset_(), columnIds_(),
243  nWeightsPerRow_(nWeightsPerRow), rowWeights_(), numNzWeight_(),
244  mayHaveDiagonalEntries(true)
245 {
246  typedef StridedData<lno_t,scalar_t> input_t;
247 
248  rowMap_ = matrix_->getRowMap();
249  colMap_ = matrix_->getColMap();
250 
251  size_t nrows = matrix_->getNodeNumRows();
252  size_t nnz = matrix_->getNodeNumEntries();
253  size_t maxnumentries =
254  matrix_->getNodeMaxNumRowEntries(); // Diff from CrsMatrix
255 
256  offset_.resize(nrows+1, 0);
257  columnIds_.resize(nnz);
258  values_.resize(nnz);
259  typename User::nonconst_local_inds_host_view_type indices("indices", maxnumentries);
260  typename User::nonconst_values_host_view_type nzs("nzs", maxnumentries);
261 
262  lno_t next = 0;
263  for (size_t i=0; i < nrows; i++){
264  lno_t row = i;
265  matrix_->getLocalRowCopy(row, indices, nzs, nnz); // Diff from CrsMatrix
266  for (size_t j=0; j < nnz; j++){
267  values_[next] = nzs[j];
268  // TODO - this will be slow
269  // Is it possible that global columns ids might be stored in order?
270  columnIds_[next++] = colMap_->getGlobalElement(indices[j]);
271  }
272  offset_[i+1] = offset_[i] + nnz;
273  }
274 
275  if (nWeightsPerRow_ > 0){
276  rowWeights_ = arcp(new input_t [nWeightsPerRow_], 0, nWeightsPerRow_, true);
277  numNzWeight_ = arcp(new bool [nWeightsPerRow_], 0, nWeightsPerRow_, true);
278  for (int i=0; i < nWeightsPerRow_; i++)
279  numNzWeight_[i] = false;
280  }
281 }
282 
284 template <typename User, typename UserCoord>
286  const scalar_t *weightVal, int stride, int idx)
287 {
288  if (this->getPrimaryEntityType() == MATRIX_ROW)
289  setRowWeights(weightVal, stride, idx);
290  else {
291  // TODO: Need to allow weights for columns and/or nonzeros
292  std::ostringstream emsg;
293  emsg << __FILE__ << "," << __LINE__
294  << " error: setWeights not yet supported for"
295  << " columns or nonzeros."
296  << std::endl;
297  throw std::runtime_error(emsg.str());
298  }
299 }
300 
302 template <typename User, typename UserCoord>
304  const scalar_t *weightVal, int stride, int idx)
305 {
306  typedef StridedData<lno_t,scalar_t> input_t;
307  if(idx<0 || idx >= nWeightsPerRow_)
308  {
309  std::ostringstream emsg;
310  emsg << __FILE__ << ":" << __LINE__
311  << " Invalid row weight index " << idx << std::endl;
312  throw std::runtime_error(emsg.str());
313  }
314 
315  size_t nvtx = getLocalNumRows();
316  ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride, false);
317  rowWeights_[idx] = input_t(weightV, stride);
318 }
319 
321 template <typename User, typename UserCoord>
323  int idx)
324 {
325  if (this->getPrimaryEntityType() == MATRIX_ROW)
326  setRowWeightIsNumberOfNonZeros(idx);
327  else {
328  // TODO: Need to allow weights for columns and/or nonzeros
329  std::ostringstream emsg;
330  emsg << __FILE__ << "," << __LINE__
331  << " error: setWeightIsNumberOfNonZeros not yet supported for"
332  << " columns" << std::endl;
333  throw std::runtime_error(emsg.str());
334  }
335 }
336 
338 template <typename User, typename UserCoord>
340  int idx)
341 {
342  if(idx<0 || idx >= nWeightsPerRow_)
343  {
344  std::ostringstream emsg;
345  emsg << __FILE__ << ":" << __LINE__
346  << " Invalid row weight index " << idx << std::endl;
347  throw std::runtime_error(emsg.str());
348  }
349 
350 
351  numNzWeight_[idx] = true;
352 }
353 
355 template <typename User, typename UserCoord>
356  template <typename Adapter>
358  const User &in, User *&out,
359  const PartitioningSolution<Adapter> &solution) const
360 {
361  // Get an import list (rows to be received)
362  size_t numNewRows;
363  ArrayRCP<gno_t> importList;
364  try{
365  numNewRows = Zoltan2::getImportList<Adapter,
367  (solution, this, importList);
368  }
370 
371  // Move the rows, creating a new matrix.
372  RCP<User> outPtr = doMigration(in, numNewRows, importList.getRawPtr());
373  out = outPtr.get();
374  outPtr.release();
375 }
376 
378 template <typename User, typename UserCoord>
379  template <typename Adapter>
381  const User &in, RCP<User> &out,
382  const PartitioningSolution<Adapter> &solution) const
383 {
384  // Get an import list (rows to be received)
385  size_t numNewRows;
386  ArrayRCP<gno_t> importList;
387  try{
388  numNewRows = Zoltan2::getImportList<Adapter,
390  (solution, this, importList);
391  }
393 
394  // Move the rows, creating a new matrix.
395  out = doMigration(in, numNewRows, importList.getRawPtr());
396 }
397 
398 
400 template < typename User, typename UserCoord>
402  const User &from,
403  size_t numLocalRows,
404  const gno_t *myNewRows
405 ) const
406 {
407  typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
408  typedef Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> tcrsmatrix_t;
409 
410  // We cannot create a Tpetra::RowMatrix, unless the underlying type is
411  // something we know (like Tpetra::CrsMatrix).
412  // If the underlying type is something different, the user probably doesn't
413  // want a Tpetra::CrsMatrix back, so we throw an error.
414 
415  // Try to cast "from" matrix to a TPetra::CrsMatrix
416  // If that fails we throw an error.
417  // We could cast as a ref which will throw std::bad_cast but with ptr
418  // approach it might be clearer what's going on here
419  const tcrsmatrix_t *pCrsMatrix = dynamic_cast<const tcrsmatrix_t *>(&from);
420 
421  if(!pCrsMatrix) {
422  throw std::logic_error("TpetraRowMatrixAdapter cannot migrate data for "
423  "your RowMatrix; it can migrate data only for "
424  "Tpetra::CrsMatrix. "
425  "You can inherit from TpetraRowMatrixAdapter and "
426  "implement migration for your RowMatrix.");
427  }
428 
429  // source map
430  const RCP<const map_t> &smap = from.getRowMap();
431  gno_t numGlobalRows = smap->getGlobalNumElements();
432  gno_t base = smap->getMinAllGlobalIndex();
433 
434  // target map
435  ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
436  const RCP<const Teuchos::Comm<int> > &comm = from.getComm();
437  RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
438 
439  // importer
440  Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
441 
442  // target matrix
443  // Chris Siefert proposed using the following to make migration
444  // more efficient.
445  // By default, the Domain and Range maps are the same as in "from".
446  // As in the original code, we instead set them both to tmap.
447  // The assumption is a square matrix.
448  // TODO: what about rectangular matrices?
449  // TODO: Should choice of domain/range maps be an option to this function?
450 
451  // KDD 3/7/16: disabling Chris' new code to avoid dashboard failures;
452  // KDD 3/7/16: can re-enable when issue #114 is fixed.
453  // KDD 3/7/16: when re-enable CSIEFERT code, can comment out
454  // KDD 3/7/16: "Original way" code.
455  // CSIEFERT RCP<tcrsmatrix_t> M;
456  // CSIEFERT from.importAndFillComplete(M, importer, tmap, tmap);
457 
458  // Original way we did it:
459  //
460  int oldNumElts = smap->getNodeNumElements();
461  int newNumElts = numLocalRows;
462 
463  // number of non zeros in my new rows
464  typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
465  vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
466  vector_t numNew(tmap); // but ETI does not yet support that.
467  for (int lid=0; lid < oldNumElts; lid++){
468  numOld.replaceGlobalValue(smap->getGlobalElement(lid),
469  scalar_t(from.getNumEntriesInLocalRow(lid)));
470  }
471  numNew.doImport(numOld, importer, Tpetra::INSERT);
472 
473  // TODO Could skip this copy if could declare vector with scalar=size_t.
474  ArrayRCP<size_t> nnz(newNumElts);
475  if (newNumElts > 0){
476  ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
477  for (int lid=0; lid < newNumElts; lid++){
478  nnz[lid] = static_cast<size_t>(ptr[lid]);
479  }
480  }
481 
482  RCP<tcrsmatrix_t> M =
483  rcp(new tcrsmatrix_t(tmap, nnz(), Tpetra::StaticProfile));
484 
485  M->doImport(from, importer, Tpetra::INSERT);
486  M->fillComplete();
487 
488  // End of original way we did it.
489  return Teuchos::rcp_dynamic_cast<User>(M);
490 }
491 
492 } //namespace Zoltan2
493 
494 #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 MatrixAdapter interface.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
InputTraits< User >::node_t node_t
InputTraits< User >::offset_t offset_t
InputTraits< User >::part_t part_t
InputTraits< User >::scalar_t scalar_t
InputTraits< User >::lno_t lno_t
InputTraits< User >::gno_t gno_t
MatrixAdapter defines the adapter interface for matrices.
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
Provides access for Zoltan2 to Tpetra::RowMatrix data.
bool CRSViewAvailable() const
Indicates whether the MatrixAdapter implements a view of the matrix in compressed sparse row (CRS) fo...
size_t getLocalNumEntries() const
Returns the number of nonzeros on this process.
size_t getLocalNumColumns() const
Returns the number of columns on this process.
size_t getLocalNumRows() const
Returns the number of rows on this process.
TpetraRowMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds, ArrayRCP< const scalar_t > &values) const
void getRowWeightsView(const scalar_t *&weights, int &stride, int idx=0) const
Provide a pointer to the row weights, if any.
void setRowWeightIsNumberOfNonZeros(int idx)
Specify an index for which the row weight should be the global number of nonzeros in the row.
void setWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each entity of the primaryEntityType.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
void getRowIDsView(const gno_t *&rowIds) const
bool useNumNonzerosAsRowWeight(int idx) const
Indicate whether row weight with index idx should be the global number of nonzeros in the row.
int getNumWeightsPerRow() const
Returns the number of weights per row (0 or greater). Row weights may be used when partitioning matri...
void setRowWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each row.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:17
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:18
Tpetra::Map map_t
Definition: mapRemotes.cpp:16
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_offset_t offset_t
The data type to represent offsets.
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.