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