Zoltan2
Loading...
Searching...
No Matches
Zoltan2_TpetraRowMatrixAdapter.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_TPETRAROWMATRIXADAPTER_HPP_
15#define _ZOLTAN2_TPETRAROWMATRIXADAPTER_HPP_
16
20
21#include <Tpetra_RowMatrix.hpp>
22
23#include <vector>
24
25namespace Zoltan2 {
26
28
41template <typename User, typename UserCoord = User>
42class TpetraRowMatrixAdapter : public MatrixAdapter<User, UserCoord> {
43public:
44
45
46#ifndef DOXYGEN_SHOULD_SKIP_THIS
49 using lno_t = typename InputTraits<User>::lno_t;
50 using gno_t = typename InputTraits<User>::gno_t;
51 using part_t = typename InputTraits<User>::part_t;
52 using node_t = typename InputTraits<User>::node_t;
53 using device_t = typename node_t::device_type;
54 using host_t = typename Kokkos::HostSpace::memory_space;
55 using user_t = User;
56 using userCoord_t = UserCoord;
57
59#endif
60
66 TpetraRowMatrixAdapter(const RCP<const User> &inmatrix,
67 int nWeightsPerRow = 0);
68
81 void setWeights(const scalar_t *weightVal, int stride, int idx = 0);
82
92 void setWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx);
93
103 void setWeightsHost(typename Base::ConstWeightsHostView1D val, int idx);
104
120 void setRowWeights(const scalar_t *weightVal, int stride, int idx = 0);
121
134 void setRowWeightsDevice(typename Base::ConstWeightsDeviceView1D val,
135 int idx);
136
149 void setRowWeightsHost(typename Base::ConstWeightsHostView1D val,
150 int idx);
151
152
158 void setWeightIsDegree(int idx);
159
166
168// The MatrixAdapter Interface
170
171 size_t getLocalNumRows() const;
172
173 size_t getLocalNumColumns() const;
174
175 size_t getLocalNumEntries() const;
176
177 bool CRSViewAvailable() const;
178
179 void getRowIDsView(const gno_t *&rowIds) const override;
180
182 typename Base::ConstIdsHostView &rowIds) const override;
183
185 typename Base::ConstIdsDeviceView &rowIds) const override;
186
187 void getCRSView(ArrayRCP<const offset_t> &offsets,
188 ArrayRCP<const gno_t> &colIds) const;
189
191 typename Base::ConstOffsetsHostView &offsets,
192 typename Base::ConstIdsHostView &colIds) const override;
193
195 typename Base::ConstOffsetsDeviceView &offsets,
196 typename Base::ConstIdsDeviceView &colIds) const override;
197
198 void getCRSView(ArrayRCP<const offset_t> &offsets,
199 ArrayRCP<const gno_t> &colIds,
200 ArrayRCP<const scalar_t> &values) const;
201
203 typename Base::ConstOffsetsHostView &offsets,
204 typename Base::ConstIdsHostView &colIds,
205 typename Base::ConstScalarsHostView &values) const override;
206
208 typename Base::ConstOffsetsDeviceView &offsets,
209 typename Base::ConstIdsDeviceView &colIds,
210 typename Base::ConstScalarsDeviceView &values) const override;
211
213
214 void getRowWeightsView(const scalar_t *&weights, int &stride,
215 int idx = 0) const;
216
217 void getRowWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
218 int idx = 0) const;
219
221 typename Base::WeightsDeviceView &weights) const override;
222
223 void getRowWeightsHostView(typename Base::WeightsHostView1D &weights,
224 int idx = 0) const;
225
227 typename Base::WeightsHostView &weights) const override;
228
229 bool useNumNonzerosAsRowWeight(int idx) const;
230
231 template <typename Adapter>
233 const User &in, User *&out,
234 const PartitioningSolution<Adapter> &solution) const;
235
236 template <typename Adapter>
238 const User &in, RCP<User> &out,
239 const PartitioningSolution<Adapter> &solution) const;
240
241protected:
242 // Used by TpetraCrsMatrixAdapter
243 TpetraRowMatrixAdapter(int nWeightsPerRow,
244 const RCP<const User> &inmatrix)
245 : matrix_(inmatrix), nWeightsPerRow_(nWeightsPerRow) {}
246
247 RCP<const User> matrix_;
248
249 ArrayRCP<offset_t> offset_;
250 ArrayRCP<gno_t> columnIds_;
251 ArrayRCP<scalar_t> values_;
252
253 typename Base::ConstOffsetsHostView offsHost_;
254 typename Base::ConstIdsHostView colIdsHost_;
255 typename Base::ScalarsHostView valuesHost_;
256
257 typename Base::ConstOffsetsDeviceView offsDevice_;
258 typename Base::ConstIdsDeviceView colIdsDevice_;
259 typename Base::ScalarsDeviceView valuesDevice_;
260
262 ArrayRCP<StridedData<lno_t, scalar_t>> rowWeights_;
263 typename Base::WeightsDeviceView rowWeightsDevice_;
264 Kokkos::View<bool *, host_t> numNzWeight_;
265
267
268 virtual RCP<User> doMigration(const User &from, size_t numLocalRows,
269 const gno_t *myNewRows) const;
270};
271
273// Definitions
275
276template <typename User, typename UserCoord>
278 const RCP<const User> &inmatrix, int nWeightsPerRow):
279 matrix_(inmatrix), offset_(), columnIds_(),
280 nWeightsPerRow_(nWeightsPerRow), rowWeights_(),
281 mayHaveDiagonalEntries(true) {
282 using strided_t = StridedData<lno_t, scalar_t>;
283 using localInds_t = typename User::nonconst_local_inds_host_view_type;
284 using localVals_t = typename User::nonconst_values_host_view_type;
285
286 const auto nrows = matrix_->getLocalNumRows();
287 const auto nnz = matrix_->getLocalNumEntries();
288 auto maxNumEntries = matrix_->getLocalMaxNumRowEntries();
289
290 // Unfortunately we have to copy the offsets, column Ids, and vals
291 // because column Ids are not usually stored in row id order.
292
293 colIdsHost_ = typename Base::ConstIdsHostView("colIdsHost_", nnz);
294 offsHost_ = typename Base::ConstOffsetsHostView("offsHost_", nrows + 1);
295 valuesHost_ = typename Base::ScalarsHostView("valuesHost_", nnz);
296
297 localInds_t localColInds("localColInds", maxNumEntries);
298 localVals_t localVals("localVals", maxNumEntries);
299
300 for (size_t r = 0; r < nrows; r++) {
301 size_t numEntries = 0;
302 matrix_->getLocalRowCopy(r, localColInds, localVals, numEntries); // Diff from CrsGraph
303
304 offsHost_(r + 1) = offsHost_(r) + numEntries;
305 for (offset_t e = offsHost_(r), i = 0; e < offsHost_(r + 1); e++) {
306 colIdsHost_(e) = matrix_->getColMap()->getGlobalElement(localColInds(i++));
307 }
308 for (size_t j = 0; j < numEntries; j++) {
309 valuesHost_(r) = localVals[j];
310 }
311 }
312 offsDevice_ = Kokkos::create_mirror_view_and_copy(
313 typename Base::device_t(), offsHost_);
314 colIdsDevice_ = Kokkos::create_mirror_view_and_copy(
315 typename Base::device_t(), colIdsHost_);
316 valuesDevice_ = Kokkos::create_mirror_view_and_copy(
317 typename Base::device_t(), valuesHost_);
318
319 if (nWeightsPerRow_ > 0) {
321 arcp(new strided_t[nWeightsPerRow_], 0, nWeightsPerRow_, true);
322
323 rowWeightsDevice_ = typename Base::WeightsDeviceView(
324 "rowWeightsDevice_", nrows, nWeightsPerRow_);
325
326 numNzWeight_ = Kokkos::View<bool *, host_t>(
327 "numNzWeight_", nWeightsPerRow_);
328
329 for (int i = 0; i < nWeightsPerRow_; ++i) {
330 numNzWeight_(i) = false;
331 }
332 }
333}
334
336template <typename User, typename UserCoord>
338 const scalar_t *weightVal, int stride, int idx) {
339 if (this->getPrimaryEntityType() == MATRIX_ROW)
340 setRowWeights(weightVal, stride, idx);
341 else {
342 // TODO: Need to allow weights for columns and/or nonzeros
343 std::ostringstream emsg;
344 emsg << __FILE__ << "," << __LINE__
345 << " error: setWeights not yet supported for"
346 << " columns or nonzeros." << std::endl;
347 throw std::runtime_error(emsg.str());
348 }
349}
350
352template <typename User, typename UserCoord>
354 typename Base::ConstWeightsDeviceView1D val, int idx) {
355 if (this->getPrimaryEntityType() == MATRIX_ROW)
356 setRowWeightsDevice(val, idx);
357 else {
358 // TODO: Need to allow weights for columns and/or nonzeros
359 std::ostringstream emsg;
360 emsg << __FILE__ << "," << __LINE__
361 << " error: setWeights not yet supported for"
362 << " columns or nonzeros." << std::endl;
363 throw std::runtime_error(emsg.str());
364 }
365}
366
368template <typename User, typename UserCoord>
370 typename Base::ConstWeightsHostView1D val, int idx) {
371 if (this->getPrimaryEntityType() == MATRIX_ROW)
372 setRowWeightsHost(val, idx);
373 else {
374 // TODO: Need to allow weights for columns and/or nonzeros
375 std::ostringstream emsg;
376 emsg << __FILE__ << "," << __LINE__
377 << " error: setWeights not yet supported for"
378 << " columns or nonzeros." << std::endl;
379 throw std::runtime_error(emsg.str());
380 }
381}
382
384template <typename User, typename UserCoord>
386 const scalar_t *weightVal, int stride, int idx) {
387 typedef StridedData<lno_t, scalar_t> input_t;
388 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
389 "Invalid row weight index: " + std::to_string(idx));
390
391 size_t nrows = getLocalNumRows();
392 ArrayRCP<const scalar_t> weightV(weightVal, 0, nrows * stride, false);
393 rowWeights_[idx] = input_t(weightV, stride);
394}
395
397template <typename User, typename UserCoord>
399 typename Base::ConstWeightsDeviceView1D weights, int idx) {
400
401 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
402 "Invalid row weight index: " + std::to_string(idx));
403
404 Kokkos::parallel_for(
405 rowWeightsDevice_.extent(0), KOKKOS_CLASS_LAMBDA(const int rowID) {
406 rowWeightsDevice_(rowID, idx) = weights(rowID);
407 });
408
409 Kokkos::fence();
410}
411
413template <typename User, typename UserCoord>
415 typename Base::ConstWeightsHostView1D weightsHost, int idx) {
416 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
417 "Invalid row weight index: " + std::to_string(idx));
418
419 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
420 typename Base::device_t(), weightsHost);
421
422 setRowWeightsDevice(weightsDevice, idx);
423}
424
426template <typename User, typename UserCoord>
428 if (this->getPrimaryEntityType() == MATRIX_ROW)
429 setRowWeightIsNumberOfNonZeros(idx);
430 else {
431 // TODO: Need to allow weights for columns and/or nonzeros
432 std::ostringstream emsg;
433 emsg << __FILE__ << "," << __LINE__
434 << " error: setWeightIsNumberOfNonZeros not yet supported for"
435 << " columns" << std::endl;
436 throw std::runtime_error(emsg.str());
437 }
438}
439
441template <typename User, typename UserCoord>
443 int idx) {
444 if (idx < 0 || idx >= nWeightsPerRow_) {
445 std::ostringstream emsg;
446 emsg << __FILE__ << ":" << __LINE__ << " Invalid row weight index " << idx
447 << std::endl;
448 throw std::runtime_error(emsg.str());
449 }
450
451 numNzWeight_(idx) = true;
452}
453
455template <typename User, typename UserCoord>
457 return matrix_->getLocalNumRows();
458}
459
461template <typename User, typename UserCoord>
463 return matrix_->getLocalNumCols();
464}
465
467template <typename User, typename UserCoord>
469 return matrix_->getLocalNumEntries();
470}
471
473template <typename User, typename UserCoord>
475
477template <typename User, typename UserCoord>
479 ArrayView<const gno_t> rowView = matrix_->getRowMap()->getLocalElementList();
480 rowIds = rowView.getRawPtr();
481}
482
484template <typename User, typename UserCoord>
486 typename Base::ConstIdsHostView &rowIds) const {
487 auto idsDevice = matrix_->getRowMap()->getMyGlobalIndices();
488 auto tmpIds = typename Base::IdsHostView("", idsDevice.extent(0));
489
490 Kokkos::deep_copy(tmpIds, idsDevice);
491
492 rowIds = tmpIds;
493}
494
496template <typename User, typename UserCoord>
498 typename Base::ConstIdsDeviceView &rowIds) const {
499
500 auto idsDevice = matrix_->getRowMap()->getMyGlobalIndices();
501 auto tmpIds = typename Base::IdsDeviceView("", idsDevice.extent(0));
502
503 Kokkos::deep_copy(tmpIds, idsDevice);
504
505 rowIds = tmpIds;
506}
507
509template <typename User, typename UserCoord>
510void TpetraRowMatrixAdapter<User, UserCoord>::getCRSView(ArrayRCP<const offset_t> &offsets,
511 ArrayRCP<const gno_t> &colIds) const {
512 offsets = offset_;
513 colIds = columnIds_;
514}
515
517template <typename User, typename UserCoord>
519 typename Base::ConstOffsetsHostView &offsets,
520 typename Base::ConstIdsHostView &colIds) const {
521 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
522 Kokkos::deep_copy(hostOffsets, offsDevice_);
523 offsets = hostOffsets;
524
525 auto hostColIds = Kokkos::create_mirror_view(colIdsDevice_);
526 Kokkos::deep_copy(hostColIds, colIdsDevice_);
527 colIds = hostColIds;
528}
529
531template <typename User, typename UserCoord>
533 typename Base::ConstOffsetsDeviceView &offsets,
534 typename Base::ConstIdsDeviceView &colIds) const {
535 offsets = offsDevice_;
536 colIds = colIdsDevice_;
537}
538
540template <typename User, typename UserCoord>
541void TpetraRowMatrixAdapter<User, UserCoord>::getCRSView(ArrayRCP<const offset_t> &offsets,
542 ArrayRCP<const gno_t> &colIds,
543 ArrayRCP<const scalar_t> &values) const {
544 offsets = offset_;
545 colIds = columnIds_;
546 values = values_;
547}
548
550template <typename User, typename UserCoord>
552 typename Base::ConstOffsetsHostView &offsets,
553 typename Base::ConstIdsHostView &colIds,
554 typename Base::ConstScalarsHostView &values) const {
555 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
556 Kokkos::deep_copy(hostOffsets, offsDevice_);
557 offsets = hostOffsets;
558
559 auto hostColIds = Kokkos::create_mirror_view(colIdsDevice_);
560 Kokkos::deep_copy(hostColIds, colIdsDevice_);
561 colIds = hostColIds;
562
563 auto hostValues = Kokkos::create_mirror_view(valuesDevice_);
564 Kokkos::deep_copy(hostValues, valuesDevice_);
565 values = hostValues;
566}
567
569template <typename User, typename UserCoord>
571 typename Base::ConstOffsetsDeviceView &offsets,
572 typename Base::ConstIdsDeviceView &colIds,
573 typename Base::ConstScalarsDeviceView &values) const {
574 offsets = offsDevice_;
575 colIds = colIdsDevice_;
576 values = valuesDevice_;
577}
578
580template <typename User, typename UserCoord>
582
584template <typename User, typename UserCoord>
586 int idx) const {
587 if (idx < 0 || idx >= nWeightsPerRow_) {
588 std::ostringstream emsg;
589 emsg << __FILE__ << ":" << __LINE__ << " Invalid row weight index "
590 << idx << std::endl;
591 throw std::runtime_error(emsg.str());
592 }
593
594 size_t length;
595 rowWeights_[idx].getStridedList(length, weights, stride);
596}
597
599template <typename User, typename UserCoord>
601 typename Base::WeightsDeviceView1D &weights, int idx) const {
602 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
603 "Invalid row weight index.");
604
605 const auto size = rowWeightsDevice_.extent(0);
606 weights = typename Base::WeightsDeviceView1D("weights", size);
607
608 Kokkos::parallel_for(
609 size, KOKKOS_CLASS_LAMBDA(const int id) {
610 weights(id) = rowWeightsDevice_(id, idx);
611 });
612
613 Kokkos::fence();
614}
615
617template <typename User, typename UserCoord>
619 typename Base::WeightsDeviceView &weights) const {
620
621 weights = rowWeightsDevice_;
622}
623
625template <typename User, typename UserCoord>
627 typename Base::WeightsHostView1D &weights, int idx) const {
628 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
629 "Invalid row weight index.");
630
631 auto weightsDevice = typename Base::WeightsDeviceView1D(
632 "weights", rowWeightsDevice_.extent(0));
633 getRowWeightsDeviceView(weightsDevice, idx);
634
635 weights = Kokkos::create_mirror_view(weightsDevice);
636 Kokkos::deep_copy(weights, weightsDevice);
637}
638
640template <typename User, typename UserCoord>
642 typename Base::WeightsHostView &weights) const {
643
644 weights = Kokkos::create_mirror_view(rowWeightsDevice_);
645 Kokkos::deep_copy(weights, rowWeightsDevice_);
646}
647
649template <typename User, typename UserCoord>
650bool TpetraRowMatrixAdapter<User, UserCoord>::useNumNonzerosAsRowWeight(int idx) const { return numNzWeight_[idx]; }
651
653template <typename User, typename UserCoord>
654template <typename Adapter>
656 const User &in, User *&out,
657 const PartitioningSolution<Adapter> &solution) const {
658 // Get an import list (rows to be received)
659 size_t numNewRows;
660 ArrayRCP<gno_t> importList;
661 try {
662 numNewRows =
663 Zoltan2::getImportList<Adapter, TpetraRowMatrixAdapter<User, UserCoord>>(
664 solution, this, importList);
665 }
667
668 // Move the rows, creating a new matrix.
669 RCP<User> outPtr = doMigration(in, numNewRows, importList.getRawPtr());
670 out = outPtr.get();
671 outPtr.release();
672}
673
675template <typename User, typename UserCoord>
676template <typename Adapter>
678 const User &in, RCP<User> &out,
679 const PartitioningSolution<Adapter> &solution) const {
680 // Get an import list (rows to be received)
681 size_t numNewRows;
682 ArrayRCP<gno_t> importList;
683 try {
684 numNewRows =
685 Zoltan2::getImportList<Adapter, TpetraRowMatrixAdapter<User, UserCoord>>(
686 solution, this, importList);
687 }
689
690 // Move the rows, creating a new matrix.
691 out = doMigration(in, numNewRows, importList.getRawPtr());
692}
693
695template <typename User, typename UserCoord>
697 const User &from, size_t numLocalRows, const gno_t *myNewRows) const {
698 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
699 typedef Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> tcrsmatrix_t;
700
701 // We cannot create a Tpetra::RowMatrix, unless the underlying type is
702 // something we know (like Tpetra::CrsMatrix).
703 // If the underlying type is something different, the user probably doesn't
704 // want a Tpetra::CrsMatrix back, so we throw an error.
705
706 // Try to cast "from" matrix to a TPetra::CrsMatrix
707 // If that fails we throw an error.
708 // We could cast as a ref which will throw std::bad_cast but with ptr
709 // approach it might be clearer what's going on here
710 const tcrsmatrix_t *pCrsMatrix = dynamic_cast<const tcrsmatrix_t *>(&from);
711
712 if (!pCrsMatrix) {
713 throw std::logic_error("TpetraRowMatrixAdapter cannot migrate data for "
714 "your RowMatrix; it can migrate data only for "
715 "Tpetra::CrsMatrix. "
716 "You can inherit from TpetraRowMatrixAdapter and "
717 "implement migration for your RowMatrix.");
718 }
719
720 // source map
721 const RCP<const map_t> &smap = from.getRowMap();
722 gno_t numGlobalRows = smap->getGlobalNumElements();
723 gno_t base = smap->getMinAllGlobalIndex();
724
725 // target map
726 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
727 const RCP<const Teuchos::Comm<int>> &comm = from.getComm();
728 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
729
730 // importer
731 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
732
733 int oldNumElts = smap->getLocalNumElements();
734 int newNumElts = numLocalRows;
735
736 // number of non zeros in my new rows
737 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
738 vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
739 vector_t numNew(tmap); // but ETI does not yet support that.
740 for (int lid = 0; lid < oldNumElts; lid++) {
741 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
742 scalar_t(from.getNumEntriesInLocalRow(lid)));
743 }
744 numNew.doImport(numOld, importer, Tpetra::INSERT);
745
746 // TODO Could skip this copy if could declare vector with scalar=size_t.
747 ArrayRCP<size_t> nnz(newNumElts);
748 if (newNumElts > 0) {
749 ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
750 for (int lid = 0; lid < newNumElts; lid++) {
751 nnz[lid] = static_cast<size_t>(ptr[lid]);
752 }
753 }
754
755 RCP<tcrsmatrix_t> M = rcp(new tcrsmatrix_t(tmap, nnz()));
756
757 M->doImport(from, importer, Tpetra::INSERT);
758 M->fillComplete();
759
760 return Teuchos::rcp_dynamic_cast<User>(M);
761}
762
763} // namespace Zoltan2
764
765#endif // _ZOLTAN2_TPETRAROWMATRIXADAPTER_HPP_
#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.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::scalar_t scalar_t
typename InputTraits< User >::gno_t gno_t
typename Kokkos::HostSpace::memory_space host_t
typename InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_t
typename node_t::device_type device_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.
void getRowWeightsHostView(typename Base::WeightsHostView &weights) const override
size_t getLocalNumRows() const
Returns the number of rows on this process.
TpetraRowMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
void getRowIDsHostView(typename Base::ConstIdsHostView &rowIds) const override
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds, ArrayRCP< const scalar_t > &values) const
ArrayRCP< StridedData< lno_t, scalar_t > > rowWeights_
void getRowWeightsView(const scalar_t *&weights, int &stride, int idx=0) const
Provide a pointer to the row weights, if any.
TpetraRowMatrixAdapter(int nWeightsPerRow, const RCP< const User > &inmatrix)
void setRowWeightIsNumberOfNonZeros(int idx)
Specify an index for which the row weight should be the global number of nonzeros in the row.
void getRowWeightsHostView(typename Base::WeightsHostView1D &weights, int idx=0) const
void setWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each entity of the primaryEntityType.
void getRowIDsView(const gno_t *&rowIds) const override
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const
void setRowWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx)
Provide a device view to row weights.
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
void getCRSHostView(typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &colIds) const override
void getRowWeightsDeviceView(typename Base::WeightsDeviceView &weights) const override
void getCRSDeviceView(typename Base::ConstOffsetsDeviceView &offsets, typename Base::ConstIdsDeviceView &colIds) const override
bool useNumNonzerosAsRowWeight(int idx) const
Indicate whether row weight with index idx should be the global number of nonzeros in the row.
void setRowWeightsHost(typename Base::ConstWeightsHostView1D val, int idx)
Provide a host view to row weights.
void setWeightsHost(typename Base::ConstWeightsHostView1D val, int idx)
Provide a host view of weights for the primary entity type.
void getCRSDeviceView(typename Base::ConstOffsetsDeviceView &offsets, typename Base::ConstIdsDeviceView &colIds, typename Base::ConstScalarsDeviceView &values) const override
int getNumWeightsPerRow() const
Returns the number of weights per row (0 or greater). Row weights may be used when partitioning matri...
void getRowWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int idx=0) const
void getCRSHostView(typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &colIds, typename Base::ConstScalarsHostView &values) const override
virtual RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows) const
void setRowWeights(const scalar_t *weightVal, int stride, int idx=0)
Specify a weight for each row.
void getRowIDsDeviceView(typename Base::ConstIdsDeviceView &rowIds) const override
void applyPartitioningSolution(const User &in, RCP< User > &out, const PartitioningSolution< Adapter > &solution) const
void setWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx)
Provide a device view of weights for the primary entity type.
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
map_t::local_ordinal_type lno_t
map_t::global_ordinal_type gno_t
Tpetra::Map map_t
Created by mbenlioglu on Aug 31, 2020.
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)
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.