14#ifndef _ZOLTAN2_TPETRAROWGRAPHADAPTER_HPP_
15#define _ZOLTAN2_TPETRAROWGRAPHADAPTER_HPP_
17#include "Kokkos_DualView.hpp"
18#include "Kokkos_UnorderedMap.hpp"
19#include <Tpetra_RowGraph.hpp>
47template <
typename User,
typename UserCoord = User>
51#ifndef DOXYGEN_SHOULD_SKIP_THIS
59 using userCoord_t = UserCoord;
74 int nEdgeWeights = 0);
226 const gno_t *&adjIds)
const override;
230 typename Base::ConstIdsDeviceView &adjIds)
const override;
233 typename Base::ConstIdsHostView &adjIds)
const override;
238 int idx)
const override;
241 int idx = 0)
const override;
244 typename Base::WeightsDeviceView &
weights)
const override;
247 int idx = 0)
const override;
250 typename Base::WeightsHostView &
weights)
const override;
257 int idx)
const override;
260 int idx = 0)
const override;
263 typename Base::WeightsDeviceView &
weights)
const override;
266 int idx = 0)
const override;
269 typename Base::WeightsHostView &
weights)
const override;
271 template <
typename Adapter>
273 const User &in, User *&out,
276 template <
typename Adapter>
278 const User &in, RCP<User> &out,
284 const RCP<const User> &graph)
305 virtual RCP<User>
doMigration(
const User &from,
size_t numLocalRows,
306 const gno_t *myNewRows)
const;
313template <
typename User,
typename UserCoord>
315 const RCP<const User> &ingraph,
int nVtxWgts,
int nEdgeWgts)
316 : graph_(ingraph), nWeightsPerVertex_(nVtxWgts),
317 nWeightsPerEdge_(nEdgeWgts), edgeWeights_() {
319 using localInds_t =
typename User::nonconst_local_inds_host_view_type;
321 const auto nvtx =
graph_->getLocalNumRows();
322 const auto nedges =
graph_->getLocalNumEntries();
324 const auto maxNumEntries =
graph_->getLocalMaxNumRowEntries();
329 adjIdsHost_ =
typename Base::ConstIdsHostView(
"adjIdsHost_", nedges);
330 offsHost_ =
typename Base::ConstOffsetsHostView(
"offsHost_", nvtx + 1);
332 localInds_t nbors(
"nbors", maxNumEntries);
334 for (
size_t v = 0; v < nvtx; v++) {
335 size_t numColInds = 0;
336 graph_->getLocalRowCopy(v, nbors, numColInds);
347 Kokkos::create_mirror_view_and_copy(
typename Base::device_t(),
offsHost_);
348 adjIdsDevice_ = Kokkos::create_mirror_view_and_copy(
typename Base::device_t(),
376template <
typename User,
typename UserCoord>
378 const scalar_t *weightVal,
int stride,
int idx) {
380 setVertexWeights(weightVal, stride, idx);
382 setEdgeWeights(weightVal, stride, idx);
386template <
typename User,
typename UserCoord>
388 typename Base::ConstWeightsDeviceView1D val,
int idx) {
390 setVertexWeightsDevice(val, idx);
392 setEdgeWeightsDevice(val, idx);
396template <
typename User,
typename UserCoord>
398 typename Base::ConstWeightsHostView1D val,
int idx) {
400 setVertexWeightsHost(val, idx);
402 setEdgeWeightsHost(val, idx);
406template <
typename User,
typename UserCoord>
408 const scalar_t *weightVal,
int stride,
int idx) {
410 "Invalid vertex weight index: " + std::to_string(idx));
412 size_t nvtx = getLocalNumVertices();
413 ArrayRCP<const scalar_t> weightV(weightVal, 0, nvtx * stride,
false);
414 vertexWeights_[idx] = input_t(weightV, stride);
418template <
typename User,
typename UserCoord>
420 typename Base::ConstWeightsDeviceView1D
weights,
int idx) {
423 "Invalid vertex weight index: " + std::to_string(idx));
428 auto vertexWeightsDevice = this->vertexWeightsDevice_;
429 Kokkos::parallel_for(
430 vertexWeightsDevice.extent(0), KOKKOS_LAMBDA(
const int vertexID) {
431 vertexWeightsDevice(vertexID, idx) =
weights(vertexID);
438template <
typename User,
typename UserCoord>
440 typename Base::ConstWeightsHostView1D weightsHost,
int idx) {
442 "Invalid vertex weight index: " + std::to_string(idx));
444 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
445 typename Base::device_t(), weightsHost);
447 setVertexWeightsDevice(weightsDevice, idx);
451template <
typename User,
typename UserCoord>
454 "setWeightIsNumberOfNonZeros is supported only for vertices");
456 setVertexWeightIsDegree(idx);
460template <
typename User,
typename UserCoord>
463 "Invalid vertex weight index.");
465 vertexDegreeWeightsHost_(idx) =
true;
469template <
typename User,
typename UserCoord>
471 const scalar_t *weightVal,
int stride,
int idx) {
475 "Invalid edge weight index" + std::to_string(idx));
477 size_t nedges = getLocalNumEdges();
478 ArrayRCP<const scalar_t> weightV(weightVal, 0, nedges * stride,
false);
479 edgeWeights_[idx] = input_t(weightV, stride);
483template <
typename User,
typename UserCoord>
485 typename Base::ConstWeightsDeviceView1D
weights,
int idx) {
487 "Invalid edge weight index.");
492 auto edgeWeightsDevice = this->edgeWeightsDevice_;
493 Kokkos::parallel_for(
494 edgeWeightsDevice.extent(0), KOKKOS_LAMBDA(
const int vertexID) {
495 edgeWeightsDevice(vertexID, idx) =
weights(vertexID);
502template <
typename User,
typename UserCoord>
504 typename Base::ConstWeightsHostView1D weightsHost,
int idx) {
506 "Invalid edge weight index.");
508 auto weightsDevice = Kokkos::create_mirror_view_and_copy(
509 typename Base::device_t(), weightsHost);
511 setEdgeWeightsDevice(weightsDevice);
515template <
typename User,
typename UserCoord>
517 return graph_->getLocalNumRows();
521template <
typename User,
typename UserCoord>
523 const gno_t *&ids)
const {
525 if (getLocalNumVertices())
526 ids = graph_->getRowMap()->getLocalElementList().getRawPtr();
530template <
typename User,
typename UserCoord>
532 typename Base::ConstIdsDeviceView &ids)
const {
536 auto idsDevice = graph_->getRowMap()->getMyGlobalIndices();
537 auto tmpIds =
typename Base::IdsDeviceView(
"", idsDevice.extent(0));
539 Kokkos::deep_copy(tmpIds, idsDevice);
545template <
typename User,
typename UserCoord>
547 typename Base::ConstIdsHostView &ids)
const {
550 auto idsDevice = graph_->getRowMap()->getMyGlobalIndices();
551 auto tmpIds =
typename Base::IdsHostView(
"", idsDevice.extent(0));
553 Kokkos::deep_copy(tmpIds, idsDevice);
559template <
typename User,
typename UserCoord>
561 return graph_->getLocalNumEntries();
565template <
typename User,
typename UserCoord>
568 offsets = offsHost_.data();
569 adjIds = (getLocalNumEdges() ? adjIdsHost_.data() : NULL);
573template <
typename User,
typename UserCoord>
575 typename Base::ConstOffsetsDeviceView &offsets,
576 typename Base::ConstIdsDeviceView &adjIds)
const {
578 offsets = offsDevice_;
579 adjIds = adjIdsDevice_;
583template <
typename User,
typename UserCoord>
585 typename Base::ConstOffsetsHostView &offsets,
586 typename Base::ConstIdsHostView &adjIds)
const {
588 auto hostIDs = Kokkos::create_mirror_view(adjIdsDevice_);
589 Kokkos::deep_copy(hostIDs, adjIdsDevice_);
592 auto hostOffsets = Kokkos::create_mirror_view(offsDevice_);
593 Kokkos::deep_copy(hostOffsets, offsDevice_);
594 offsets = hostOffsets;
598template <
typename User,
typename UserCoord>
600 return nWeightsPerVertex_;
604template <
typename User,
typename UserCoord>
609 "Invalid vertex weight index.");
612 vertexWeights_[idx].getStridedList(length,
weights, stride);
616template <
typename User,
typename UserCoord>
618 typename Base::WeightsDeviceView1D &
weights,
int idx)
const {
620 "Invalid vertex weight index.");
622 const auto size = vertexWeightsDevice_.extent(0);
623 weights =
typename Base::WeightsDeviceView1D(
"weights", size);
625 auto vertexWeightsDevice = this->vertexWeightsDevice_;
626 Kokkos::parallel_for(
627 size, KOKKOS_LAMBDA(
const int id) {
628 weights(
id) = vertexWeightsDevice(
id, idx);
635template <
typename User,
typename UserCoord>
637 typename Base::WeightsDeviceView &
weights)
const {
639 weights = vertexWeightsDevice_;
643template <
typename User,
typename UserCoord>
645 typename Base::WeightsHostView1D &
weights,
int idx)
const {
647 "Invalid vertex weight index.");
649 auto weightsDevice =
typename Base::WeightsDeviceView1D(
650 "weights", vertexWeightsDevice_.extent(0));
651 getVertexWeightsDeviceView(weightsDevice, idx);
653 weights = Kokkos::create_mirror_view(weightsDevice);
654 Kokkos::deep_copy(
weights, weightsDevice);
658template <
typename User,
typename UserCoord>
660 typename Base::WeightsHostView &
weights)
const {
662 weights = Kokkos::create_mirror_view(vertexWeightsDevice_);
663 Kokkos::deep_copy(
weights, vertexWeightsDevice_);
667template <
typename User,
typename UserCoord>
670 return vertexDegreeWeightsHost_(idx);
674template <
typename User,
typename UserCoord>
676 return nWeightsPerEdge_;
680template <
typename User,
typename UserCoord>
684 "Invalid edge weight index.");
687 edgeWeights_[idx].getStridedList(length,
weights, stride);
691template <
typename User,
typename UserCoord>
693 typename Base::WeightsDeviceView1D &
weights,
int idx)
const {
695 weights = Kokkos::subview(edgeWeightsDevice_, Kokkos::ALL, idx);
699template <
typename User,
typename UserCoord>
701 typename Base::WeightsDeviceView &
weights)
const {
707template <
typename User,
typename UserCoord>
709 typename Base::WeightsHostView1D &
weights,
int idx)
const {
711 auto weightsDevice = Kokkos::subview(edgeWeightsDevice_, Kokkos::ALL, idx);
712 weights = Kokkos::create_mirror_view(weightsDevice);
713 Kokkos::deep_copy(
weights, weightsDevice);
717template <
typename User,
typename UserCoord>
719 typename Base::WeightsHostView &
weights)
const {
721 weights = Kokkos::create_mirror_view(edgeWeightsDevice_);
722 Kokkos::deep_copy(
weights, edgeWeightsDevice_);
726template <
typename User,
typename UserCoord>
727template <
typename Adapter>
729 const User &in, User *&out,
733 ArrayRCP<gno_t> importList;
736 Zoltan2::getImportList<Adapter, TpetraRowGraphAdapter<User, UserCoord>>(
737 solution,
this, importList);
742 RCP<User> outPtr = doMigration(in, numNewVtx, importList.getRawPtr());
748template <
typename User,
typename UserCoord>
749template <
typename Adapter>
751 const User &in, RCP<User> &out,
755 ArrayRCP<gno_t> importList;
758 Zoltan2::getImportList<Adapter, TpetraRowGraphAdapter<User, UserCoord>>(
759 solution,
this, importList);
764 out = doMigration(in, numNewVtx, importList.getRawPtr());
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;
783 const tcrsgraph_t *pCrsGraphSrc =
dynamic_cast<const tcrsgraph_t *
>(&from);
786 throw std::logic_error(
"TpetraRowGraphAdapter cannot migrate data for "
787 "your RowGraph; it can migrate data only for "
789 "You can inherit from TpetraRowGraphAdapter and "
790 "implement migration for your RowGraph.");
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();
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));
805 Tpetra::Import<lno_t, gno_t, node_t> importer(smap, tmap);
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));
815 numNew.doImport(numOld, importer, Tpetra::INSERT);
817 size_t numElts = tmap->getLocalNumElements();
818 ArrayRCP<const gno_t> nnz;
820 nnz = numNew.getData(0);
822 ArrayRCP<const size_t> nnz_size_t;
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]);
831 nnz_size_t = arcp_reinterpret_cast<const size_t>(nnz);
835 RCP<tcrsgraph_t> G = rcp(
new tcrsgraph_t(tmap, nnz_size_t()));
837 G->doImport(*pCrsGraphSrc, importer, Tpetra::INSERT);
839 return Teuchos::rcp_dynamic_cast<User>(G);
#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.
Base::WeightsDeviceView edgeWeightsDevice_
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.
Base::ConstOffsetsDeviceView offsDevice_
Base::WeightsDeviceView vertexWeightsDevice_
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
Base::ConstIdsDeviceView adjIdsDevice_
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
Base::ConstOffsetsHostView offsHost_
Base::ConstIdsHostView adjIdsHost_
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.
Base::VtxDegreeHostView vertexDegreeWeightsHost_
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
Created by mbenlioglu on Aug 31, 2020.
static void AssertCondition(bool condition, const std::string &message, const char *file=__FILE__, int line=__LINE__)