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