Zoltan2
Loading...
Searching...
No Matches
Zoltan2_XpetraCrsMatrixAdapter.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Zoltan2: A package of combinatorial algorithms for scientific computing
4//
5// Copyright 2012 NTESS and the Zoltan2 contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14#ifndef _ZOLTAN2_XPETRACRSMATRIXADAPTER_HPP_
15#define _ZOLTAN2_XPETRACRSMATRIXADAPTER_HPP_
16
21
22#include <Xpetra_CrsMatrix.hpp>
23
24#include <iostream>
25#include <cassert>
26
27namespace Zoltan2 {
28
30
51template <typename User, typename UserCoord=User>
52 class XpetraCrsMatrixAdapter : public MatrixAdapter<User, UserCoord> {
53public:
54
55#ifndef DOXYGEN_SHOULD_SKIP_THIS
57 using lno_t = typename InputTraits<User>::lno_t;
58 using gno_t = typename InputTraits<User>::gno_t;
59 using part_t = typename InputTraits<User>::part_t;
60 using node_t = typename InputTraits<User>::node_t;
62 using xmatrix_t = Xpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t>;
63
64 using userCoord_t = UserCoord;
65 using user_t = User;
66#endif
67
71
77 XpetraCrsMatrixAdapter(const RCP<const User> &inmatrix,
78 int nWeightsPerRow=0);
79
92 void setWeights(const scalar_t *weightVal, int stride, int idx = 0);
93
109 void setRowWeights(const scalar_t *weightVal, int stride, int idx = 0);
110
116 void setWeightIsDegree(int idx);
117
123 void setRowWeightIsNumberOfNonZeros(int idx);
124
126 // The MatrixAdapter interface.
128
129 size_t getLocalNumRows() const {
130 return matrix_->getLocalNumRows();
131 }
132
133 size_t getLocalNumColumns() const {
134 return matrix_->getLocalNumCols();
135 }
136
137 size_t getLocalNumEntries() const {
138 return matrix_->getLocalNumEntries();
139 }
140
141 void getRowIDsView(const gno_t *&rowIds) const
142 {
143 ArrayView<const gno_t> rowView = rowMap_->getLocalElementList();
144 rowIds = rowView.getRawPtr();
145 }
146
147 void getColumnIDsView(const gno_t *&colIds) const
148 {
149 ArrayView<const gno_t> colView = colMap_->getLocalElementList();
150 colIds = colView.getRawPtr();
151 }
152
153 void getCRSView(ArrayRCP<const offset_t> &offsets, ArrayRCP<const gno_t> &colIds) const
154 {
155 ArrayRCP< const lno_t > localColumnIds;
156 ArrayRCP<const scalar_t> values;
157 matrix_->getAllValues(offsets,localColumnIds,values);
158 colIds = columnIds_;
159 }
160
161 bool CRSViewAvailable() const { return true; }
162
163 void getCRSView(ArrayRCP<const offset_t> &offsets,
164 ArrayRCP<const gno_t> &colIds,
165 ArrayRCP<const scalar_t> &values) const {
166 ArrayRCP< const lno_t > localColumnIds;
167 matrix_->getAllValues(offsets,localColumnIds,values);
168 colIds = columnIds_;
169 }
170
171 void getCCSView(ArrayRCP<const offset_t> &offsets,
172 ArrayRCP<const gno_t> &rowIds) const override {
173 ArrayRCP<const offset_t> crsOffsets;
174 ArrayRCP<const lno_t> crsLocalColumnIds;
175 ArrayRCP<const scalar_t> values;
176 matrix_->getAllValues(crsOffsets, crsLocalColumnIds, values);
177
178 const auto localRowIds = rowMap_->getLocalElementList();
179 const auto numLocalCols = colMap_->getLocalNumElements();
180
181 // Lambda used to compute local row based on column index from CRS view
182 auto determineRow = [&crsOffsets, &localRowIds](const int columnIdx) {
183 int curLocalRow = 0;
184 for (int rowIdx = 0; rowIdx < localRowIds.size(); ++rowIdx) {
185 if (rowIdx < (localRowIds.size() - 1)) {
186 if (static_cast<offset_t>(columnIdx) < crsOffsets[rowIdx + 1]) {
187 return curLocalRow;
188 }
189 ++curLocalRow;
190 } else {
191 return curLocalRow;
192 }
193 }
194
195 return -1;
196 };
197
198 // Vector of global rows per each local column
199 std::vector<std::vector<gno_t>> rowIDsPerCol(numLocalCols);
200
201 for (int colIdx = 0; colIdx < crsLocalColumnIds.size(); ++colIdx) {
202 const auto colID = crsLocalColumnIds[colIdx];
203 const auto globalRow = rowMap_->getGlobalElement(determineRow(colIdx));
204
205 rowIDsPerCol[colID].push_back(globalRow);
206 }
207
208 size_t offsetWrite = 0;
209 ArrayRCP<gno_t> ccsRowIds(values.size());
210 ArrayRCP<offset_t> ccsOffsets(colMap_->getLocalNumElements() + 1);
211
212 ccsOffsets[0] = 0;
213 for (int64_t colID = 1; colID < ccsOffsets.size(); ++colID) {
214 const auto &rowIDs = rowIDsPerCol[colID - 1];
215
216 if (not rowIDs.empty()) {
217 std::copy(rowIDs.begin(), rowIDs.end(),
218 ccsRowIds.begin() + offsetWrite);
219 offsetWrite += rowIDs.size();
220 }
221
222 ccsOffsets[colID] = offsetWrite;
223 }
224
225 ccsOffsets[numLocalCols] = crsLocalColumnIds.size();
226
227 rowIds = ccsRowIds;
228 offsets = ccsOffsets;
229 }
230
231 int getNumWeightsPerRow() const { return nWeightsPerRow_; }
232
233 void getRowWeightsView(const scalar_t *&weights, int &stride,
234 int idx = 0) const
235 {
236 if(idx<0 || idx >= nWeightsPerRow_)
237 {
238 std::ostringstream emsg;
239 emsg << __FILE__ << ":" << __LINE__
240 << " Invalid row weight index " << idx << std::endl;
241 throw std::runtime_error(emsg.str());
242 }
243
244 size_t length;
245 rowWeights_[idx].getStridedList(length, weights, stride);
246 }
247
248 bool useNumNonzerosAsRowWeight(int idx) const { return numNzWeight_[idx];}
249
250 template <typename Adapter>
251 void applyPartitioningSolution(const User &in, User *&out,
252 const PartitioningSolution<Adapter> &solution) const;
253
254 template <typename Adapter>
255 void applyPartitioningSolution(const User &in, RCP<User> &out,
256 const PartitioningSolution<Adapter> &solution) const;
257
258private:
259
260 RCP<const User> inmatrix_;
261 RCP<const xmatrix_t> matrix_;
262 RCP<const Xpetra::Map<lno_t, gno_t, node_t> > rowMap_;
263 RCP<const Xpetra::Map<lno_t, gno_t, node_t> > colMap_;
264 lno_t base_;
265 ArrayRCP<gno_t> columnIds_; // TODO: Refactor adapter to localColumnIds_
266
267 int nWeightsPerRow_;
268 ArrayRCP<StridedData<lno_t, scalar_t> > rowWeights_;
269 ArrayRCP<bool> numNzWeight_;
270
271 bool mayHaveDiagonalEntries;
272};
273
275// Definitions
277
278template <typename User, typename UserCoord>
280 const RCP<const User> &inmatrix, int nWeightsPerRow):
281 inmatrix_(inmatrix), matrix_(), rowMap_(), colMap_(),
282 columnIds_(),
283 nWeightsPerRow_(nWeightsPerRow), rowWeights_(), numNzWeight_(),
284 mayHaveDiagonalEntries(true)
285{
286 typedef StridedData<lno_t,scalar_t> input_t;
287 try {
288 matrix_ = rcp_const_cast<const xmatrix_t>(
289 XpetraTraits<User>::convertToXpetra(rcp_const_cast<User>(inmatrix)));
290 }
292
293 rowMap_ = matrix_->getRowMap();
294 colMap_ = matrix_->getColMap();
295
296 size_t nrows = matrix_->getLocalNumRows();
297 size_t nnz = matrix_->getLocalNumEntries();
298
299 // Get ArrayRCP pointers to the structures in the underlying matrix
300 ArrayRCP< const offset_t > offset;
301 ArrayRCP< const lno_t > localColumnIds;
302 ArrayRCP< const scalar_t > values;
303 matrix_->getAllValues(offset,localColumnIds,values);
304 columnIds_.resize(nnz, 0);
305
306 for (offset_t i = 0; i < offset[nrows]; i++) {
307 columnIds_[i] = colMap_->getGlobalElement(localColumnIds[i]);
308 }
309
310 if (nWeightsPerRow_ > 0){
311 rowWeights_ = arcp(new input_t [nWeightsPerRow_], 0, nWeightsPerRow_, true);
312 numNzWeight_ = arcp(new bool [nWeightsPerRow_], 0, nWeightsPerRow_, true);
313 for (int i=0; i < nWeightsPerRow_; i++)
314 numNzWeight_[i] = false;
315 }
316}
317
319template <typename User, typename UserCoord>
321 const scalar_t *weightVal, int stride, int idx)
322{
323 if (this->getPrimaryEntityType() == MATRIX_ROW)
324 setRowWeights(weightVal, stride, idx);
325 else {
326 // TODO: Need to allow weights for columns and/or nonzeros
327 std::ostringstream emsg;
328 emsg << __FILE__ << "," << __LINE__
329 << " error: setWeights not yet supported for"
330 << " columns or nonzeros."
331 << std::endl;
332 throw std::runtime_error(emsg.str());
333 }
334}
335
337template <typename User, typename UserCoord>
339 const scalar_t *weightVal, int stride, int idx)
340{
341 typedef StridedData<lno_t,scalar_t> input_t;
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 size_t nvtx = getLocalNumRows();
351 ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx*stride, false);
352 rowWeights_[idx] = input_t(weightV, stride);
353}
354
356template <typename User, typename UserCoord>
358 int idx)
359{
360 if (this->getPrimaryEntityType() == MATRIX_ROW)
361 setRowWeightIsNumberOfNonZeros(idx);
362 else {
363 // TODO: Need to allow weights for columns and/or nonzeros
364 std::ostringstream emsg;
365 emsg << __FILE__ << "," << __LINE__
366 << " error: setWeightIsNumberOfNonZeros not yet supported for"
367 << " columns" << std::endl;
368 throw std::runtime_error(emsg.str());
369 }
370}
371
373template <typename User, typename UserCoord>
375 int idx)
376{
377 if(idx<0 || idx >= nWeightsPerRow_)
378 {
379 std::ostringstream emsg;
380 emsg << __FILE__ << ":" << __LINE__
381 << " Invalid row weight index " << idx << std::endl;
382 throw std::runtime_error(emsg.str());
383 }
384
385
386 numNzWeight_[idx] = true;
387}
388
390template <typename User, typename UserCoord>
391 template <typename Adapter>
393 const User &in, User *&out,
394 const PartitioningSolution<Adapter> &solution) const
395{
396 // Get an import list (rows to be received)
397 size_t numNewRows;
398 ArrayRCP<gno_t> importList;
399 try{
400 numNewRows = Zoltan2::getImportList<Adapter,
402 (solution, this, importList);
403 }
405
406 // Move the rows, creating a new matrix.
407 RCP<User> outPtr = XpetraTraits<User>::doMigration(in, numNewRows,
408 importList.getRawPtr());
409 out = const_cast<User *>(outPtr.get());
410 outPtr.release();
411}
412
414template <typename User, typename UserCoord>
415 template <typename Adapter>
417 const User &in, RCP<User> &out,
418 const PartitioningSolution<Adapter> &solution) const
419{
420 // Get an import list (rows to be received)
421 size_t numNewRows;
422 ArrayRCP<gno_t> importList;
423 try{
424 numNewRows = Zoltan2::getImportList<Adapter,
426 (solution, this, importList);
427 }
429
430 // Move the rows, creating a new matrix.
431 out = XpetraTraits<User>::doMigration(in, numNewRows,
432 importList.getRawPtr());
433}
434
435} //namespace Zoltan2
436
437#endif
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
Xpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > xmatrix_t
#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.
Traits of Xpetra classes, including migration method.
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_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 Xpetra::CrsMatrix data.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
bool CRSViewAvailable() const
Indicates whether the MatrixAdapter implements a view of the matrix in compressed sparse row (CRS) fo...
bool useNumNonzerosAsRowWeight(int idx) const
Indicate whether row weight with index idx should be the global number of nonzeros in the row.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const
void getRowIDsView(const gno_t *&rowIds) const
size_t getLocalNumColumns() const
Returns the number of columns on this process.
void getCCSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &rowIds) const override
void getColumnIDsView(const gno_t *&colIds) const
size_t getLocalNumEntries() const
Returns the number of nonzeros on this process.
void setRowWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each row.
size_t getLocalNumRows() const
Returns the number of rows on this process.
int getNumWeightsPerRow() const
Returns the number of weights per row (0 or greater). Row weights may be used when partitioning matri...
void setWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each entity of the primaryEntityType.
void getRowWeightsView(const scalar_t *&weights, int &stride, int idx=0) const
Provide a pointer to the row weights, if any.
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds, ArrayRCP< const scalar_t > &values) const
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
void setRowWeightIsNumberOfNonZeros(int idx)
Specify an index for which the row weight should be the global number of nonzeros in the row.
XpetraCrsMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
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.
Defines the traits required for Tpetra, Eptra and Xpetra objects.
static RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows)
Migrate the object Given a user object and a new row distribution, create and return a new user objec...