50 #ifndef _ZOLTAN2_TPETRAROWMATRIXADAPTER_HPP_ 51 #define _ZOLTAN2_TPETRAROWMATRIXADAPTER_HPP_ 57 #include <Tpetra_RowMatrix.hpp> 75 template <
typename User,
typename UserCoord=User>
79 #ifndef DOXYGEN_SHOULD_SKIP_THIS 86 typedef UserCoord userCoord_t;
99 int nWeightsPerRow=0);
113 void setWeights(
const scalar_t *weightVal,
int stride,
int idx = 0);
130 void setRowWeights(
const scalar_t *weightVal,
int stride,
int idx = 0);
151 return matrix_->getNodeNumRows();
155 return matrix_->getNodeNumCols();
159 return matrix_->getNodeNumEntries();
166 ArrayView<const gno_t> rowView = rowMap_->getNodeElementList();
167 rowIds = rowView.getRawPtr();
170 void getCRSView(
const lno_t *&offsets,
const gno_t *&colIds)
const 172 offsets = offset_.getRawPtr();
173 colIds = columnIds_.getRawPtr();
177 const scalar_t *&values)
const 179 offsets = offset_.getRawPtr();
180 colIds = columnIds_.getRawPtr();
181 values = values_.getRawPtr();
190 if(idx<0 || idx >= nWeightsPerRow_)
192 std::ostringstream emsg;
193 emsg << __FILE__ <<
":" << __LINE__
194 <<
" Invalid row weight index " << idx << std::endl;
195 throw std::runtime_error(emsg.str());
200 rowWeights_[idx].getStridedList(length, weights, stride);
205 template <
typename Adapter>
209 template <
typename Adapter>
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_;
220 ArrayRCP<scalar_t> values_;
223 ArrayRCP<StridedData<lno_t, scalar_t> > rowWeights_;
224 ArrayRCP<bool> numNzWeight_;
226 bool mayHaveDiagonalEntries;
228 RCP<User> doMigration(
const User &from,
size_t numLocalRows,
229 const gno_t *myNewRows)
const;
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)
246 rowMap_ = matrix_->getRowMap();
247 colMap_ = matrix_->getColMap();
249 size_t nrows = matrix_->getNodeNumRows();
250 size_t nnz = matrix_->getNodeNumEntries();
251 size_t maxnumentries =
252 matrix_->getNodeMaxNumRowEntries();
254 offset_.resize(nrows+1, 0);
255 columnIds_.resize(nnz);
257 ArrayRCP<lno_t> indices(maxnumentries);
258 ArrayRCP<scalar_t> nzs(maxnumentries);
260 for (
size_t i=0; i < nrows; i++){
262 matrix_->getLocalRowCopy(row, indices(), nzs(), nnz);
263 for (
size_t j=0; j < nnz; j++){
264 values_[next] = nzs[j];
267 columnIds_[next++] = colMap_->getGlobalElement(indices[j]);
269 offset_[i+1] = offset_[i] + nnz;
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;
281 template <
typename User,
typename UserCoord>
283 const scalar_t *weightVal,
int stride,
int idx)
289 std::ostringstream emsg;
290 emsg << __FILE__ <<
"," << __LINE__
291 <<
" error: setWeights not yet supported for" 292 <<
" columns or nonzeros." 294 throw std::runtime_error(emsg.str());
299 template <
typename User,
typename UserCoord>
301 const scalar_t *weightVal,
int stride,
int idx)
304 if(idx<0 || idx >= nWeightsPerRow_)
306 std::ostringstream emsg;
307 emsg << __FILE__ <<
":" << __LINE__
308 <<
" Invalid row weight index " << idx << std::endl;
309 throw std::runtime_error(emsg.str());
313 ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride,
false);
314 rowWeights_[idx] = input_t(weightV, stride);
318 template <
typename User,
typename UserCoord>
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());
335 template <
typename User,
typename UserCoord>
339 if(idx<0 || idx >= nWeightsPerRow_)
341 std::ostringstream emsg;
342 emsg << __FILE__ <<
":" << __LINE__
343 <<
" Invalid row weight index " << idx << std::endl;
344 throw std::runtime_error(emsg.str());
348 numNzWeight_[idx] =
true;
352 template <
typename User,
typename UserCoord>
353 template <
typename Adapter>
355 const User &in, User *&out,
360 ArrayRCP<gno_t> importList;
364 (solution,
this, importList);
369 RCP<User> outPtr = doMigration(in, numNewRows, importList.getRawPtr());
375 template <
typename User,
typename UserCoord>
376 template <
typename Adapter>
378 const User &in, RCP<User> &out,
383 ArrayRCP<gno_t> importList;
387 (solution,
this, importList);
392 out = doMigration(in, numNewRows, importList.getRawPtr());
397 template <
typename User,
typename UserCoord>
401 const gno_t *myNewRows
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;
416 const tcrsmatrix_t *pCrsMatrix =
dynamic_cast<const tcrsmatrix_t *
>(&from);
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.");
427 const RCP<const map_t> &smap = from.getRowMap();
428 gno_t numGlobalRows = smap->getGlobalNumElements();
429 gno_t base = smap->getMinAllGlobalIndex();
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));
437 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
457 int oldNumElts = smap->getNodeNumElements();
458 int newNumElts = numLocalRows;
461 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
462 vector_t numOld(smap);
463 vector_t numNew(tmap);
464 for (
int lid=0; lid < oldNumElts; lid++){
465 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
466 scalar_t(from.getNumEntriesInLocalRow(lid)));
468 numNew.doImport(numOld, importer, Tpetra::INSERT);
471 ArrayRCP<size_t> nnz(newNumElts);
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]);
479 RCP<tcrsmatrix_t> M = rcp(
new tcrsmatrix_t(tmap, nnz,
480 Tpetra::StaticProfile));
481 M->doImport(from, importer, Tpetra::INSERT);
485 return Teuchos::rcp_dynamic_cast<User>(M);
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...
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...
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...
~TpetraRowMatrixAdapter()
Destructor.
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.
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.
void getCRSView(const lno_t *&offsets, const gno_t *&colIds) const
Sets pointers to this process' 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...
void getRowIDsView(const gno_t *&rowIds) const
Sets pointer to this process' rows' global IDs.
TpetraRowMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
This file defines the StridedData class.