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 override;
172
173 size_t getLocalNumColumns() const override;
174
175 size_t getLocalNumEntries() const override;
176
177 bool CRSViewAvailable() const override;
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 override;
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 override;
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
212 int getNumWeightsPerRow() const override;
213
214 void getRowWeightsView(const scalar_t *&weights, int &stride,
215 int idx = 0) const override;
216
217 void getRowWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
218 int idx = 0) const override;
219
221 typename Base::WeightsDeviceView &weights) const override;
222
223 void getRowWeightsHostView(typename Base::WeightsHostView1D &weights,
224 int idx = 0) const override;
225
227 typename Base::WeightsHostView &weights) const override;
228
229 bool useNumNonzerosAsRowWeight(int idx) const override;
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 auto rowWeightsDevice = this->rowWeightsDevice_;
405 Kokkos::parallel_for(
406 rowWeightsDevice.extent(0), KOKKOS_LAMBDA(const int rowID) {
407 rowWeightsDevice(rowID, idx) = weights(rowID);
408 });
409
410 Kokkos::fence();
411}
412
414template <typename User, typename UserCoord>
416 typename Base::ConstWeightsHostView1D weightsHost, int idx) {
417 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
418 "Invalid row weight index: " + std::to_string(idx));
419
420 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
421 typename Base::device_t(), weightsHost);
422
423 setRowWeightsDevice(weightsDevice, idx);
424}
425
427template <typename User, typename UserCoord>
429 if (this->getPrimaryEntityType() == MATRIX_ROW)
430 setRowWeightIsNumberOfNonZeros(idx);
431 else {
432 // TODO: Need to allow weights for columns and/or nonzeros
433 std::ostringstream emsg;
434 emsg << __FILE__ << "," << __LINE__
435 << " error: setWeightIsNumberOfNonZeros not yet supported for"
436 << " columns" << std::endl;
437 throw std::runtime_error(emsg.str());
438 }
439}
440
442template <typename User, typename UserCoord>
444 int idx) {
445 if (idx < 0 || idx >= nWeightsPerRow_) {
446 std::ostringstream emsg;
447 emsg << __FILE__ << ":" << __LINE__ << " Invalid row weight index " << idx
448 << std::endl;
449 throw std::runtime_error(emsg.str());
450 }
451
452 numNzWeight_(idx) = true;
453}
454
456template <typename User, typename UserCoord>
458 return matrix_->getLocalNumRows();
459}
460
462template <typename User, typename UserCoord>
464 return matrix_->getLocalNumCols();
465}
466
468template <typename User, typename UserCoord>
470 return matrix_->getLocalNumEntries();
471}
472
474template <typename User, typename UserCoord>
476
478template <typename User, typename UserCoord>
480 ArrayView<const gno_t> rowView = matrix_->getRowMap()->getLocalElementList();
481 rowIds = rowView.getRawPtr();
482}
483
485template <typename User, typename UserCoord>
487 typename Base::ConstIdsHostView &rowIds) const {
488 auto idsDevice = matrix_->getRowMap()->getMyGlobalIndices();
489 auto tmpIds = typename Base::IdsHostView("", idsDevice.extent(0));
490
491 Kokkos::deep_copy(tmpIds, idsDevice);
492
493 rowIds = tmpIds;
494}
495
497template <typename User, typename UserCoord>
499 typename Base::ConstIdsDeviceView &rowIds) const {
500
501 auto idsDevice = matrix_->getRowMap()->getMyGlobalIndices();
502 auto tmpIds = typename Base::IdsDeviceView("", idsDevice.extent(0));
503
504 Kokkos::deep_copy(tmpIds, idsDevice);
505
506 rowIds = tmpIds;
507}
508
510template <typename User, typename UserCoord>
511void TpetraRowMatrixAdapter<User, UserCoord>::getCRSView(ArrayRCP<const offset_t> &offsets,
512 ArrayRCP<const gno_t> &colIds) const {
513 offsets = offset_;
514 colIds = columnIds_;
515}
516
518template <typename User, typename UserCoord>
520 typename Base::ConstOffsetsHostView &offsets,
521 typename Base::ConstIdsHostView &colIds) const {
522 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
523 Kokkos::deep_copy(hostOffsets, offsDevice_);
524 offsets = hostOffsets;
525
526 auto hostColIds = Kokkos::create_mirror_view(colIdsDevice_);
527 Kokkos::deep_copy(hostColIds, colIdsDevice_);
528 colIds = hostColIds;
529}
530
532template <typename User, typename UserCoord>
534 typename Base::ConstOffsetsDeviceView &offsets,
535 typename Base::ConstIdsDeviceView &colIds) const {
536 offsets = offsDevice_;
537 colIds = colIdsDevice_;
538}
539
541template <typename User, typename UserCoord>
542void TpetraRowMatrixAdapter<User, UserCoord>::getCRSView(ArrayRCP<const offset_t> &offsets,
543 ArrayRCP<const gno_t> &colIds,
544 ArrayRCP<const scalar_t> &values) const {
545 offsets = offset_;
546 colIds = columnIds_;
547 values = values_;
548}
549
551template <typename User, typename UserCoord>
553 typename Base::ConstOffsetsHostView &offsets,
554 typename Base::ConstIdsHostView &colIds,
555 typename Base::ConstScalarsHostView &values) const {
556 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
557 Kokkos::deep_copy(hostOffsets, offsDevice_);
558 offsets = hostOffsets;
559
560 auto hostColIds = Kokkos::create_mirror_view(colIdsDevice_);
561 Kokkos::deep_copy(hostColIds, colIdsDevice_);
562 colIds = hostColIds;
563
564 auto hostValues = Kokkos::create_mirror_view(valuesDevice_);
565 Kokkos::deep_copy(hostValues, valuesDevice_);
566 values = hostValues;
567}
568
570template <typename User, typename UserCoord>
572 typename Base::ConstOffsetsDeviceView &offsets,
573 typename Base::ConstIdsDeviceView &colIds,
574 typename Base::ConstScalarsDeviceView &values) const {
575 offsets = offsDevice_;
576 colIds = colIdsDevice_;
577 values = valuesDevice_;
578}
579
581template <typename User, typename UserCoord>
583
585template <typename User, typename UserCoord>
587 int idx) const {
588 if (idx < 0 || idx >= nWeightsPerRow_) {
589 std::ostringstream emsg;
590 emsg << __FILE__ << ":" << __LINE__ << " Invalid row weight index "
591 << idx << std::endl;
592 throw std::runtime_error(emsg.str());
593 }
594
595 size_t length;
596 rowWeights_[idx].getStridedList(length, weights, stride);
597}
598
600template <typename User, typename UserCoord>
602 typename Base::WeightsDeviceView1D &weights, int idx) const {
603 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
604 "Invalid row weight index.");
605
606 const auto size = rowWeightsDevice_.extent(0);
607 weights = typename Base::WeightsDeviceView1D("weights", size);
608
609 auto rowWeightsDevice = this->rowWeightsDevice_;
610 Kokkos::parallel_for(
611 size, KOKKOS_LAMBDA(const int id) {
612 weights(id) = rowWeightsDevice(id, idx);
613 });
614
615 Kokkos::fence();
616}
617
619template <typename User, typename UserCoord>
621 typename Base::WeightsDeviceView &weights) const {
622
623 weights = rowWeightsDevice_;
624}
625
627template <typename User, typename UserCoord>
629 typename Base::WeightsHostView1D &weights, int idx) const {
630 AssertCondition((idx >= 0) and (idx < nWeightsPerRow_),
631 "Invalid row weight index.");
632
633 auto weightsDevice = typename Base::WeightsDeviceView1D(
634 "weights", rowWeightsDevice_.extent(0));
635 getRowWeightsDeviceView(weightsDevice, idx);
636
637 weights = Kokkos::create_mirror_view(weightsDevice);
638 Kokkos::deep_copy(weights, weightsDevice);
639}
640
642template <typename User, typename UserCoord>
644 typename Base::WeightsHostView &weights) const {
645
646 weights = Kokkos::create_mirror_view(rowWeightsDevice_);
647 Kokkos::deep_copy(weights, rowWeightsDevice_);
648}
649
651template <typename User, typename UserCoord>
652bool TpetraRowMatrixAdapter<User, UserCoord>::useNumNonzerosAsRowWeight(int idx) const { return numNzWeight_[idx]; }
653
655template <typename User, typename UserCoord>
656template <typename Adapter>
658 const User &in, User *&out,
659 const PartitioningSolution<Adapter> &solution) const {
660 // Get an import list (rows to be received)
661 size_t numNewRows;
662 ArrayRCP<gno_t> importList;
663 try {
664 numNewRows =
665 Zoltan2::getImportList<Adapter, TpetraRowMatrixAdapter<User, UserCoord>>(
666 solution, this, importList);
667 }
669
670 // Move the rows, creating a new matrix.
671 RCP<User> outPtr = doMigration(in, numNewRows, importList.getRawPtr());
672 out = outPtr.get();
673 outPtr.release();
674}
675
677template <typename User, typename UserCoord>
678template <typename Adapter>
680 const User &in, RCP<User> &out,
681 const PartitioningSolution<Adapter> &solution) const {
682 // Get an import list (rows to be received)
683 size_t numNewRows;
684 ArrayRCP<gno_t> importList;
685 try {
686 numNewRows =
687 Zoltan2::getImportList<Adapter, TpetraRowMatrixAdapter<User, UserCoord>>(
688 solution, this, importList);
689 }
691
692 // Move the rows, creating a new matrix.
693 out = doMigration(in, numNewRows, importList.getRawPtr());
694}
695
697template <typename User, typename UserCoord>
699 const User &from, size_t numLocalRows, const gno_t *myNewRows) const {
700 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
701 typedef Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t> tcrsmatrix_t;
702
703 // We cannot create a Tpetra::RowMatrix, unless the underlying type is
704 // something we know (like Tpetra::CrsMatrix).
705 // If the underlying type is something different, the user probably doesn't
706 // want a Tpetra::CrsMatrix back, so we throw an error.
707
708 // Try to cast "from" matrix to a TPetra::CrsMatrix
709 // If that fails we throw an error.
710 // We could cast as a ref which will throw std::bad_cast but with ptr
711 // approach it might be clearer what's going on here
712 const tcrsmatrix_t *pCrsMatrix = dynamic_cast<const tcrsmatrix_t *>(&from);
713
714 if (!pCrsMatrix) {
715 throw std::logic_error("TpetraRowMatrixAdapter cannot migrate data for "
716 "your RowMatrix; it can migrate data only for "
717 "Tpetra::CrsMatrix. "
718 "You can inherit from TpetraRowMatrixAdapter and "
719 "implement migration for your RowMatrix.");
720 }
721
722 // source map
723 const RCP<const map_t> &smap = from.getRowMap();
724 gno_t numGlobalRows = smap->getGlobalNumElements();
725 gno_t base = smap->getMinAllGlobalIndex();
726
727 // target map
728 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
729 const RCP<const Teuchos::Comm<int>> &comm = from.getComm();
730 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
731
732 // importer
733 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
734
735 int oldNumElts = smap->getLocalNumElements();
736 int newNumElts = numLocalRows;
737
738 // number of non zeros in my new rows
739 typedef Tpetra::Vector<scalar_t, lno_t, gno_t, node_t> vector_t;
740 vector_t numOld(smap); // TODO These vectors should have scalar=size_t,
741 vector_t numNew(tmap); // but ETI does not yet support that.
742 for (int lid = 0; lid < oldNumElts; lid++) {
743 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
744 scalar_t(from.getNumEntriesInLocalRow(lid)));
745 }
746 numNew.doImport(numOld, importer, Tpetra::INSERT);
747
748 // TODO Could skip this copy if could declare vector with scalar=size_t.
749 ArrayRCP<size_t> nnz(newNumElts);
750 if (newNumElts > 0) {
751 ArrayRCP<scalar_t> ptr = numNew.getDataNonConst(0);
752 for (int lid = 0; lid < newNumElts; lid++) {
753 nnz[lid] = static_cast<size_t>(ptr[lid]);
754 }
755 }
756
757 RCP<tcrsmatrix_t> M = rcp(new tcrsmatrix_t(tmap, nnz()));
758
759 M->doImport(from, importer, Tpetra::INSERT);
760 M->fillComplete();
761
762 return Teuchos::rcp_dynamic_cast<User>(M);
763}
764
765} // namespace Zoltan2
766
767#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.
size_t getLocalNumRows() const override
Returns the number of rows on this process.
void getRowWeightsHostView(typename Base::WeightsHostView &weights) const override
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds, ArrayRCP< const scalar_t > &values) const override
TpetraRowMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
void getRowIDsHostView(typename Base::ConstIdsHostView &rowIds) const override
size_t getLocalNumColumns() const override
Returns the number of columns on this process.
ArrayRCP< StridedData< lno_t, scalar_t > > rowWeights_
TpetraRowMatrixAdapter(int nWeightsPerRow, const RCP< const User > &inmatrix)
bool CRSViewAvailable() const override
Indicates whether the MatrixAdapter implements a view of the matrix in compressed sparse row (CRS) fo...
void setRowWeightIsNumberOfNonZeros(int idx)
Specify an index for which the row weight should be the global number of nonzeros in the row.
size_t getLocalNumEntries() const override
Returns the number of nonzeros on this process.
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 getRowWeightsHostView(typename Base::WeightsHostView1D &weights, int idx=0) const override
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.
int getNumWeightsPerRow() const override
Returns the number of weights per row (0 or greater). Row weights may be used when partitioning matri...
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
void getRowWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int idx=0) const override
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
void getCRSHostView(typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &colIds, typename Base::ConstScalarsHostView &values) const override
void getRowWeightsView(const scalar_t *&weights, int &stride, int idx=0) const override
Provide a pointer to the row weights, if any.
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
void getCRSView(ArrayRCP< const offset_t > &offsets, ArrayRCP< const gno_t > &colIds) const override
bool useNumNonzerosAsRowWeight(int idx) const override
Indicate whether row weight with index idx should be the global number of nonzeros in the row.
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.