Zoltan2
Loading...
Searching...
No Matches
Zoltan2_TpetraRowGraphAdapter.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_TPETRAROWGRAPHADAPTER_HPP_
15#define _ZOLTAN2_TPETRAROWGRAPHADAPTER_HPP_
16
17#include "Kokkos_DualView.hpp"
18#include "Kokkos_UnorderedMap.hpp"
19#include <Tpetra_RowGraph.hpp>
23#include <string>
24
25namespace Zoltan2 {
26
47template <typename User, typename UserCoord = User>
48class TpetraRowGraphAdapter : public GraphAdapter<User, UserCoord> {
49
50public:
51#ifndef DOXYGEN_SHOULD_SKIP_THIS
54 using lno_t = typename InputTraits<User>::lno_t;
55 using gno_t = typename InputTraits<User>::gno_t;
56 using part_t = typename InputTraits<User>::part_t;
57 using node_t = typename InputTraits<User>::node_t;
58 using user_t = User;
59 using userCoord_t = UserCoord;
60
62#endif
63
73 TpetraRowGraphAdapter(const RCP<const User> &ingraph, int nVtxWeights = 0,
74 int nEdgeWeights = 0);
75
88 void setWeights(const scalar_t *val, int stride, int idx);
89
98 void setWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx);
99
108 void setWeightsHost(typename Base::ConstWeightsHostView1D val, int idx);
109
125 void setVertexWeights(const scalar_t *val, int stride, int idx);
126
138 void setVertexWeightsDevice(typename Base::ConstWeightsDeviceView1D val,
139 int idx);
140
152 void setVertexWeightsHost(typename Base::ConstWeightsHostView1D val, int idx);
153
159 void setWeightIsDegree(int idx);
160
167
190 void setEdgeWeights(const scalar_t *val, int stride, int idx);
191
197 void setEdgeWeightsDevice(typename Base::ConstWeightsDeviceView1D val,
198 int idx);
199
205 void setEdgeWeightsHost(typename Base::ConstWeightsHostView1D val, int idx);
206
208 // The GraphAdapter interface.
210
211 // TODO: Assuming rows == objects;
212 // TODO: Need to add option for columns or nonzeros?
213 size_t getLocalNumVertices() const override;
214
215 void getVertexIDsView(const gno_t *&ids) const override;
216
217 void
218 getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override;
219
220 void
221 getVertexIDsHostView(typename Base::ConstIdsHostView &ids) const override;
222
223 size_t getLocalNumEdges() const override;
224
225 void getEdgesView(const offset_t *&offsets,
226 const gno_t *&adjIds) const override;
227
228 void
229 getEdgesDeviceView(typename Base::ConstOffsetsDeviceView &offsets,
230 typename Base::ConstIdsDeviceView &adjIds) const override;
231
232 void getEdgesHostView(typename Base::ConstOffsetsHostView &offsets,
233 typename Base::ConstIdsHostView &adjIds) const override;
234
235 int getNumWeightsPerVertex() const override;
236
237 void getVertexWeightsView(const scalar_t *&weights, int &stride,
238 int idx) const override;
239
240 void getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
241 int idx = 0) const override;
242
244 typename Base::WeightsDeviceView &weights) const override;
245
246 void getVertexWeightsHostView(typename Base::WeightsHostView1D &weights,
247 int idx = 0) const override;
248
250 typename Base::WeightsHostView &weights) const override;
251
252 bool useDegreeAsVertexWeight(int idx) const override;
253
254 int getNumWeightsPerEdge() const override;
255
256 void getEdgeWeightsView(const scalar_t *&weights, int &stride,
257 int idx) const override;
258
259 void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights,
260 int idx = 0) const override;
261
263 typename Base::WeightsDeviceView &weights) const override;
264
265 void getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights,
266 int idx = 0) const override;
267
269 typename Base::WeightsHostView &weights) const override;
270
271 template <typename Adapter>
273 const User &in, User *&out,
274 const PartitioningSolution<Adapter> &solution) const;
275
276 template <typename Adapter>
278 const User &in, RCP<User> &out,
279 const PartitioningSolution<Adapter> &solution) const;
280
281protected:
282 // Useb by TpetraCrsGraphAdapter
283 TpetraRowGraphAdapter(int nVtxWgts, int nEdgeWgts,
284 const RCP<const User> &graph)
285 : graph_(graph), nWeightsPerVertex_(nVtxWgts),
286 nWeightsPerEdge_(nEdgeWgts) {}
287
288 RCP<const User> graph_;
289
290 typename Base::ConstOffsetsHostView offsHost_;
291 typename Base::ConstIdsHostView adjIdsHost_;
292
293 typename Base::ConstIdsDeviceView adjIdsDevice_;
294 typename Base::ConstOffsetsDeviceView offsDevice_;
295
297 ArrayRCP<StridedData<lno_t, scalar_t>> vertexWeights_;
298 typename Base::WeightsDeviceView vertexWeightsDevice_;
299 typename Base::VtxDegreeHostView vertexDegreeWeightsHost_;
300
302 ArrayRCP<StridedData<lno_t, scalar_t>> edgeWeights_;
303 typename Base::WeightsDeviceView edgeWeightsDevice_;
304
305 virtual RCP<User> doMigration(const User &from, size_t numLocalRows,
306 const gno_t *myNewRows) const;
307};
308
310// Definitions
312
313template <typename User, typename UserCoord>
315 const RCP<const User> &ingraph, int nVtxWgts, int nEdgeWgts)
316 : graph_(ingraph), nWeightsPerVertex_(nVtxWgts),
317 nWeightsPerEdge_(nEdgeWgts), edgeWeights_() {
318 using strided_t = StridedData<lno_t, scalar_t>;
319 using localInds_t = typename User::nonconst_local_inds_host_view_type;
320
321 const auto nvtx = graph_->getLocalNumRows();
322 const auto nedges = graph_->getLocalNumEntries();
323 // Diff from CrsMatrix
324 const auto maxNumEntries = graph_->getLocalMaxNumRowEntries();
325
326 // Unfortunately we have to copy the offsets and edge Ids
327 // because edge Ids are not usually stored in vertex id order.
328
329 adjIdsHost_ = typename Base::ConstIdsHostView("adjIdsHost_", nedges);
330 offsHost_ = typename Base::ConstOffsetsHostView("offsHost_", nvtx + 1);
331
332 localInds_t nbors("nbors", maxNumEntries);
333
334 for (size_t v = 0; v < nvtx; v++) {
335 size_t numColInds = 0;
336 graph_->getLocalRowCopy(v, nbors, numColInds); // Diff from CrsGraph
337
338 offsHost_(v + 1) = offsHost_(v) + numColInds;
339 for (offset_t e = offsHost_(v), i = 0; e < offsHost_(v + 1); e++) {
340 adjIdsHost_(e) = graph_->getColMap()->getGlobalElement(nbors(i++));
341 }
342 }
343
344 // Since there's no direct getter of offsets and edges in device view,
345 // we have to deep copy here
347 Kokkos::create_mirror_view_and_copy(typename Base::device_t(), offsHost_);
348 adjIdsDevice_ = Kokkos::create_mirror_view_and_copy(typename Base::device_t(),
350
351 if (nWeightsPerVertex_ > 0) {
353 arcp(new strided_t[nWeightsPerVertex_], 0, nWeightsPerVertex_, true);
354
355 vertexWeightsDevice_ = typename Base::WeightsDeviceView(
356 "vertexWeightsDevice_", nvtx, nWeightsPerVertex_);
357
358 vertexDegreeWeightsHost_ = typename Base::VtxDegreeHostView(
359 "vertexDegreeWeightsHost_", nWeightsPerVertex_);
360
361 for (int i = 0; i < nWeightsPerVertex_; ++i) {
362 vertexDegreeWeightsHost_(i) = false;
363 }
364 }
365
366 if (nWeightsPerEdge_ > 0) {
368 arcp(new strided_t[nWeightsPerEdge_], 0, nWeightsPerEdge_, true);
369
370 edgeWeightsDevice_ = typename Base::WeightsDeviceView(
371 "nWeightsPerEdge_", graph_->getLocalNumRows(), nWeightsPerEdge_);
372 }
373}
374
376template <typename User, typename UserCoord>
378 const scalar_t *weightVal, int stride, int idx) {
379 if (this->getPrimaryEntityType() == GRAPH_VERTEX)
380 setVertexWeights(weightVal, stride, idx);
381 else
382 setEdgeWeights(weightVal, stride, idx);
383}
384
386template <typename User, typename UserCoord>
388 typename Base::ConstWeightsDeviceView1D val, int idx) {
389 if (this->getPrimaryEntityType() == GRAPH_VERTEX)
390 setVertexWeightsDevice(val, idx);
391 else
392 setEdgeWeightsDevice(val, idx);
393}
394
396template <typename User, typename UserCoord>
398 typename Base::ConstWeightsHostView1D val, int idx) {
399 if (this->getPrimaryEntityType() == GRAPH_VERTEX)
400 setVertexWeightsHost(val, idx);
401 else
402 setEdgeWeightsHost(val, idx);
403}
404
406template <typename User, typename UserCoord>
408 const scalar_t *weightVal, int stride, int idx) {
409 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
410 "Invalid vertex weight index: " + std::to_string(idx));
411
412 size_t nvtx = getLocalNumVertices();
413 ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx * stride, false);
414 vertexWeights_[idx] = input_t(weightV, stride);
415}
416
418template <typename User, typename UserCoord>
420 typename Base::ConstWeightsDeviceView1D weights, int idx) {
421
422 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
423 "Invalid vertex weight index: " + std::to_string(idx));
424
425 AssertCondition(vertexWeightsDevice_.extent(0) == weights.extent(0),
426 "Invalid sizes!");
427
428 auto vertexWeightsDevice = this->vertexWeightsDevice_;
429 Kokkos::parallel_for(
430 vertexWeightsDevice.extent(0), KOKKOS_LAMBDA(const int vertexID) {
431 vertexWeightsDevice(vertexID, idx) = weights(vertexID);
432 });
433
434 Kokkos::fence();
435}
436
438template <typename User, typename UserCoord>
440 typename Base::ConstWeightsHostView1D weightsHost, int idx) {
441 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
442 "Invalid vertex weight index: " + std::to_string(idx));
443
444 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
445 typename Base::device_t(), weightsHost);
446
447 setVertexWeightsDevice(weightsDevice, idx);
448}
449
451template <typename User, typename UserCoord>
453 AssertCondition(this->getPrimaryEntityType() == GRAPH_VERTEX,
454 "setWeightIsNumberOfNonZeros is supported only for vertices");
455
456 setVertexWeightIsDegree(idx);
457}
458
460template <typename User, typename UserCoord>
462 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
463 "Invalid vertex weight index.");
464
465 vertexDegreeWeightsHost_(idx) = true;
466}
467
469template <typename User, typename UserCoord>
471 const scalar_t *weightVal, int stride, int idx) {
472 typedef StridedData<lno_t, scalar_t> input_t;
473
474 AssertCondition((idx >= 0) and (idx < nWeightsPerEdge_),
475 "Invalid edge weight index" + std::to_string(idx));
476
477 size_t nedges = getLocalNumEdges();
478 ArrayRCP<const scalar_t> weightV(weightVal, 0, nedges * stride, false);
479 edgeWeights_[idx] = input_t(weightV, stride);
480}
481
483template <typename User, typename UserCoord>
485 typename Base::ConstWeightsDeviceView1D weights, int idx) {
486 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
487 "Invalid edge weight index.");
488
489 AssertCondition(edgeWeightsDevice_.extent(0) == weights.extent(0),
490 "Invalid sizes!");
491
492 auto edgeWeightsDevice = this->edgeWeightsDevice_;
493 Kokkos::parallel_for(
494 edgeWeightsDevice.extent(0), KOKKOS_LAMBDA(const int vertexID) {
495 edgeWeightsDevice(vertexID, idx) = weights(vertexID);
496 });
497
498 Kokkos::fence();
499}
500
502template <typename User, typename UserCoord>
504 typename Base::ConstWeightsHostView1D weightsHost, int idx) {
505 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
506 "Invalid edge weight index.");
507
508 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
509 typename Base::device_t(), weightsHost);
510
511 setEdgeWeightsDevice(weightsDevice);
512}
513
515template <typename User, typename UserCoord>
517 return graph_->getLocalNumRows();
518}
519
521template <typename User, typename UserCoord>
523 const gno_t *&ids) const {
524 ids = NULL;
525 if (getLocalNumVertices())
526 ids = graph_->getRowMap()->getLocalElementList().getRawPtr();
527}
528
530template <typename User, typename UserCoord>
532 typename Base::ConstIdsDeviceView &ids) const {
533
534 // TODO: Making a ConstIdsDeviceView LayoutLeft would proably remove the
535 // need of creating tmpIds
536 auto idsDevice = graph_->getRowMap()->getMyGlobalIndices();
537 auto tmpIds = typename Base::IdsDeviceView("", idsDevice.extent(0));
538
539 Kokkos::deep_copy(tmpIds, idsDevice);
540
541 ids = tmpIds;
542}
543
545template <typename User, typename UserCoord>
547 typename Base::ConstIdsHostView &ids) const {
548 // TODO: Making a ConstIdsDeviceView LayoutLeft would proably remove the
549 // need of creating tmpIds
550 auto idsDevice = graph_->getRowMap()->getMyGlobalIndices();
551 auto tmpIds = typename Base::IdsHostView("", idsDevice.extent(0));
552
553 Kokkos::deep_copy(tmpIds, idsDevice);
554
555 ids = tmpIds;
556}
557
559template <typename User, typename UserCoord>
561 return graph_->getLocalNumEntries();
562}
563
565template <typename User, typename UserCoord>
567 const offset_t *&offsets, const gno_t *&adjIds) const {
568 offsets = offsHost_.data();
569 adjIds = (getLocalNumEdges() ? adjIdsHost_.data() : NULL);
570}
571
573template <typename User, typename UserCoord>
575 typename Base::ConstOffsetsDeviceView &offsets,
576 typename Base::ConstIdsDeviceView &adjIds) const {
577
578 offsets = offsDevice_;
579 adjIds = adjIdsDevice_;
580}
581
583template <typename User, typename UserCoord>
585 typename Base::ConstOffsetsHostView &offsets,
586 typename Base::ConstIdsHostView &adjIds) const {
587
588 auto hostIDs = Kokkos::create_mirror_view(adjIdsDevice_);
589 Kokkos::deep_copy(hostIDs, adjIdsDevice_);
590 adjIds = hostIDs;
591
592 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
593 Kokkos::deep_copy(hostOffsets, offsDevice_);
594 offsets = hostOffsets;
595}
596
598template <typename User, typename UserCoord>
600 return nWeightsPerVertex_;
601}
602
604template <typename User, typename UserCoord>
606 const scalar_t *&weights, int &stride, int idx) const {
607
608 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
609 "Invalid vertex weight index.");
610
611 size_t length;
612 vertexWeights_[idx].getStridedList(length, weights, stride);
613}
614
616template <typename User, typename UserCoord>
618 typename Base::WeightsDeviceView1D &weights, int idx) const {
619 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
620 "Invalid vertex weight index.");
621
622 const auto size = vertexWeightsDevice_.extent(0);
623 weights = typename Base::WeightsDeviceView1D("weights", size);
624
625 auto vertexWeightsDevice = this->vertexWeightsDevice_;
626 Kokkos::parallel_for(
627 size, KOKKOS_LAMBDA(const int id) {
628 weights(id) = vertexWeightsDevice(id, idx);
629 });
630
631 Kokkos::fence();
632}
633
635template <typename User, typename UserCoord>
637 typename Base::WeightsDeviceView &weights) const {
638
639 weights = vertexWeightsDevice_;
640}
641
643template <typename User, typename UserCoord>
645 typename Base::WeightsHostView1D &weights, int idx) const {
646 AssertCondition((idx >= 0) and (idx < nWeightsPerVertex_),
647 "Invalid vertex weight index.");
648
649 auto weightsDevice = typename Base::WeightsDeviceView1D(
650 "weights", vertexWeightsDevice_.extent(0));
651 getVertexWeightsDeviceView(weightsDevice, idx);
652
653 weights = Kokkos::create_mirror_view(weightsDevice);
654 Kokkos::deep_copy(weights, weightsDevice);
655}
656
658template <typename User, typename UserCoord>
660 typename Base::WeightsHostView &weights) const {
661
662 weights = Kokkos::create_mirror_view(vertexWeightsDevice_);
663 Kokkos::deep_copy(weights, vertexWeightsDevice_);
664}
665
667template <typename User, typename UserCoord>
669 int idx) const {
670 return vertexDegreeWeightsHost_(idx);
671}
672
674template <typename User, typename UserCoord>
676 return nWeightsPerEdge_;
677}
678
680template <typename User, typename UserCoord>
682 const scalar_t *&weights, int &stride, int idx) const {
683 AssertCondition((idx >= 0) and (idx < nWeightsPerEdge_),
684 "Invalid edge weight index.");
685
686 size_t length;
687 edgeWeights_[idx].getStridedList(length, weights, stride);
688}
689
691template <typename User, typename UserCoord>
693 typename Base::WeightsDeviceView1D &weights, int idx) const {
694
695 weights = Kokkos::subview(edgeWeightsDevice_, Kokkos::ALL, idx);
696}
697
699template <typename User, typename UserCoord>
701 typename Base::WeightsDeviceView &weights) const {
702
703 weights = edgeWeightsDevice_;
704}
705
707template <typename User, typename UserCoord>
709 typename Base::WeightsHostView1D &weights, int idx) const {
710
711 auto weightsDevice = Kokkos::subview(edgeWeightsDevice_, Kokkos::ALL, idx);
712 weights = Kokkos::create_mirror_view(weightsDevice);
713 Kokkos::deep_copy(weights, weightsDevice);
714}
715
717template <typename User, typename UserCoord>
719 typename Base::WeightsHostView &weights) const {
720
721 weights = Kokkos::create_mirror_view(edgeWeightsDevice_);
722 Kokkos::deep_copy(weights, edgeWeightsDevice_);
723}
724
726template <typename User, typename UserCoord>
727template <typename Adapter>
729 const User &in, User *&out,
730 const PartitioningSolution<Adapter> &solution) const {
731 // Get an import list (rows to be received)
732 size_t numNewVtx;
733 ArrayRCP<gno_t> importList;
734 try {
735 numNewVtx =
736 Zoltan2::getImportList<Adapter, TpetraRowGraphAdapter<User, UserCoord>>(
737 solution, this, importList);
738 }
740
741 // Move the rows, creating a new graph.
742 RCP<User> outPtr = doMigration(in, numNewVtx, importList.getRawPtr());
743 out = outPtr.get();
744 outPtr.release();
745}
746
748template <typename User, typename UserCoord>
749template <typename Adapter>
751 const User &in, RCP<User> &out,
752 const PartitioningSolution<Adapter> &solution) const {
753 // Get an import list (rows to be received)
754 size_t numNewVtx;
755 ArrayRCP<gno_t> importList;
756 try {
757 numNewVtx =
758 Zoltan2::getImportList<Adapter, TpetraRowGraphAdapter<User, UserCoord>>(
759 solution, this, importList);
760 }
762
763 // Move the rows, creating a new graph.
764 out = doMigration(in, numNewVtx, importList.getRawPtr());
765}
766
768template <typename User, typename UserCoord>
770 const User &from, size_t numLocalRows, const gno_t *myNewRows) const {
771 typedef Tpetra::Map<lno_t, gno_t, node_t> map_t;
772 typedef Tpetra::CrsGraph<lno_t, gno_t, node_t> tcrsgraph_t;
773
774 // We cannot create a Tpetra::RowGraph, unless the underlying type is
775 // something we know (like Tpetra::CrsGraph).
776 // If the underlying type is something different, the user probably doesn't
777 // want a Tpetra::CrsGraph back, so we throw an error.
778
779 // Try to cast "from" graph to a TPetra::CrsGraph
780 // If that fails we throw an error.
781 // We could cast as a ref which will throw std::bad_cast but with ptr
782 // approach it might be clearer what's going on here
783 const tcrsgraph_t *pCrsGraphSrc = dynamic_cast<const tcrsgraph_t *>(&from);
784
785 if (!pCrsGraphSrc) {
786 throw std::logic_error("TpetraRowGraphAdapter cannot migrate data for "
787 "your RowGraph; it can migrate data only for "
788 "Tpetra::CrsGraph. "
789 "You can inherit from TpetraRowGraphAdapter and "
790 "implement migration for your RowGraph.");
791 }
792
793 // source map
794 const RCP<const map_t> &smap = from.getRowMap();
795 int oldNumElts = smap->getLocalNumElements();
796 gno_t numGlobalRows = smap->getGlobalNumElements();
797 gno_t base = smap->getMinAllGlobalIndex();
798
799 // target map
800 ArrayView<const gno_t> rowList(myNewRows, numLocalRows);
801 const RCP<const Teuchos::Comm<int>> &comm = from.getComm();
802 RCP<const map_t> tmap = rcp(new map_t(numGlobalRows, rowList, base, comm));
803
804 // importer
805 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
806
807 // number of entries in my new rows
808 typedef Tpetra::Vector<gno_t, lno_t, gno_t, node_t> vector_t;
809 vector_t numOld(smap);
810 vector_t numNew(tmap);
811 for (int lid = 0; lid < oldNumElts; lid++) {
812 numOld.replaceGlobalValue(smap->getGlobalElement(lid),
813 from.getNumEntriesInLocalRow(lid));
814 }
815 numNew.doImport(numOld, importer, Tpetra::INSERT);
816
817 size_t numElts = tmap->getLocalNumElements();
818 ArrayRCP<const gno_t> nnz;
819 if (numElts > 0)
820 nnz = numNew.getData(0); // hangs if vector len == 0
821
822 ArrayRCP<const size_t> nnz_size_t;
823
824 if (numElts && sizeof(gno_t) != sizeof(size_t)) {
825 size_t *vals = new size_t[numElts];
826 nnz_size_t = arcp(vals, 0, numElts, true);
827 for (size_t i = 0; i < numElts; i++) {
828 vals[i] = static_cast<size_t>(nnz[i]);
829 }
830 } else {
831 nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
832 }
833
834 // target graph
835 RCP<tcrsgraph_t> G = rcp(new tcrsgraph_t(tmap, nnz_size_t()));
836
837 G->doImport(*pCrsGraphSrc, importer, Tpetra::INSERT);
838 G->fillComplete();
839 return Teuchos::rcp_dynamic_cast<User>(G);
840}
841
842} // namespace Zoltan2
843
844#endif
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Defines the GraphAdapter 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 InputTraits< User >::offset_t offset_t
typename InputTraits< User >::part_t part_t
GraphAdapter defines the interface for graph-based user data.
A PartitioningSolution is a solution to a partitioning problem.
The StridedData class manages lists of weights or coordinates.
Provides access for Zoltan2 to Tpetra::RowGraph data.
void getEdgeWeightsHostView(typename Base::WeightsHostView &weights) const override
Provide a host view of the edge weights, if any.
ArrayRCP< StridedData< lno_t, scalar_t > > vertexWeights_
void setEdgeWeightsHost(typename Base::ConstWeightsHostView1D val, int idx)
Provide a host view to edge weights.
ArrayRCP< StridedData< lno_t, scalar_t > > edgeWeights_
void getEdgesHostView(typename Base::ConstOffsetsHostView &offsets, typename Base::ConstIdsHostView &adjIds) const override
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
size_t getLocalNumVertices() const override
Returns the number of vertices on this process.
void getEdgesView(const offset_t *&offsets, const gno_t *&adjIds) const override
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
void setEdgeWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to edge weights.
void getVertexWeightsHostView(typename Base::WeightsHostView &weights) const override
Provide a host view of the vertex weights, if any.
void getVertexWeightsView(const scalar_t *&weights, int &stride, int idx) const override
Provide a pointer to the vertex weights, if any.
void getVertexIDsHostView(typename Base::ConstIdsHostView &ids) const override
Sets pointers to this process' graph entries.
TpetraRowGraphAdapter(int nVtxWgts, int nEdgeWgts, const RCP< const User > &graph)
void getVertexWeightsDeviceView(typename Base::WeightsDeviceView &weights) const override
Provide a device view of the vertex weights, if any.
void getVertexWeightsHostView(typename Base::WeightsHostView1D &weights, int idx=0) const override
Provide a host view of the vertex weights, if any.
void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView &weights) const override
Provide a device view of the edge weights, if any.
void applyPartitioningSolution(const User &in, RCP< User > &out, const PartitioningSolution< Adapter > &solution) const
void setWeightIsDegree(int idx)
Specify an index for which the weight should be the degree of the entity.
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 getEdgeWeightsView(const scalar_t *&weights, int &stride, int idx) const override
Provide a pointer to the edge weights, if any.
void getVertexWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int idx=0) const override
Provide a device view of the vertex weights, if any.
void getVertexIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override
Sets pointers to this process' graph entries.
virtual RCP< User > doMigration(const User &from, size_t numLocalRows, const gno_t *myNewRows) const
void getEdgeWeightsHostView(typename Base::WeightsHostView1D &weights, int idx=0) const override
Provide a host view of the edge weights, if any.
void setWeightsHost(typename Base::ConstWeightsHostView1D val, int idx)
Provide a host view of weights for the primary entity type.
void setWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to weights for the primary entity type.
void getVertexIDsView(const gno_t *&ids) const override
Sets pointers to this process' graph entries.
void setVertexWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx)
Provide a device view to vertex weights.
int getNumWeightsPerEdge() const override
Returns the number (0 or greater) of edge weights.
int getNumWeightsPerVertex() const override
Returns the number (0 or greater) of weights per vertex.
size_t getLocalNumEdges() const override
Returns the number of edges on this process.
void setVertexWeights(const scalar_t *val, int stride, int idx)
Provide a pointer to vertex weights.
TpetraRowGraphAdapter(const RCP< const User > &ingraph, int nVtxWeights=0, int nEdgeWeights=0)
Constructor for graph with no weights or coordinates.
bool useDegreeAsVertexWeight(int idx) const override
Indicate whether vertex weight with index idx should be the global degree of the vertex.
void setVertexWeightIsDegree(int idx)
Specify an index for which the vertex weight should be the degree of the vertex.
void setEdgeWeightsDevice(typename Base::ConstWeightsDeviceView1D val, int idx)
Provide a device view to edge weights.
void getEdgeWeightsDeviceView(typename Base::WeightsDeviceView1D &weights, int idx=0) const override
Provide a device view of the edge weights, if any.
void getEdgesDeviceView(typename Base::ConstOffsetsDeviceView &offsets, typename Base::ConstIdsDeviceView &adjIds) const override
Gets adjacency lists for all vertices in a compressed sparse row (CSR) format.
void setVertexWeightsHost(typename Base::ConstWeightsHostView1D val, int idx)
Provide a host view to vertex weights.
map_t::local_ordinal_type lno_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.