Zoltan2
Loading...
Searching...
No Matches
Zoltan2_TpetraCrsMatrixAdapter.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_TPETRACRSMATRIXADAPTER_HPP_
15#define _ZOLTAN2_TPETRACRSMATRIXADAPTER_HPP_
16
22
23#include <Tpetra_CrsMatrix.hpp>
24
25#include <vector>
26
27namespace Zoltan2 {
28
30
49template <typename User, typename UserCoord = User>
50 class TpetraCrsMatrixAdapter : public TpetraRowMatrixAdapter<User,UserCoord> {
51
52public:
53#ifndef DOXYGEN_SHOULD_SKIP_THIS
55 using lno_t = typename InputTraits<User>::lno_t;
56 using gno_t = typename InputTraits<User>::gno_t;
57 using part_t = typename InputTraits<User>::part_t;
58 using node_t = typename InputTraits<User>::node_t;
60 using tmatrix_t = Tpetra::CrsMatrix<scalar_t, lno_t, gno_t, node_t>;
61 using device_t = typename node_t::device_type;
62 using host_t = typename Kokkos::HostSpace::memory_space;
63 using user_t = User;
64 using userCoord_t = UserCoord;
65
68 #endif
69
75 TpetraCrsMatrixAdapter(const RCP<const User> &inmatrix,
76 int nWeightsPerRow=0);
77
78 void init(const RCP<const User> &inmatrix);
79
82 RCP<const User> getUserMatrix() const { return this->matrix_; }
83
84 template <typename Adapter>
85 void applyPartitioningSolution(const User &in, User *&out,
86 const PartitioningSolution<Adapter> &solution) const;
87
88 template <typename Adapter>
89 void applyPartitioningSolution(const User &in, RCP<User> &out,
90 const PartitioningSolution<Adapter> &solution) const;
91};
92
94// Definitions
96
97 template <typename User, typename UserCoord>
98 void TpetraCrsMatrixAdapter<User,UserCoord>::init(const RCP<const User> &inmatrix) {
99 auto colIdsDevice = inmatrix->getLocalIndicesDevice();
100
101 auto colIdsGlobalDevice =
102 typename TpetraCrsMatrixAdapter<User, UserCoord>::Base::IdsDeviceView("colIdsGlobalDevice", colIdsDevice.extent(0));
103 auto colMap = inmatrix->getColMap();
104 auto lclColMap = colMap->getLocalMap();
105
106 // Convert to global IDs using Tpetra::Map
107 Kokkos::parallel_for("colIdsGlobalDevice",
108 Kokkos::RangePolicy<typename User::node_type::execution_space>(
109 0, colIdsGlobalDevice.extent(0)),
110 KOKKOS_LAMBDA(const int i) {
111 colIdsGlobalDevice(i) =
112 lclColMap.getGlobalElement(colIdsDevice(i));
113 });
114
115 this->colIdsDevice_ = colIdsGlobalDevice;
116 this->offsDevice_ = inmatrix->getLocalRowPtrsDevice();
117 }
118
119 template <typename User, typename UserCoord>
121 const RCP<const User> &inmatrix, int nWeightsPerRow):
122 RowMatrix(nWeightsPerRow, inmatrix) {
123
124 init(inmatrix);
125
126 if (this->nWeightsPerRow_ > 0) {
127
128 this->rowWeightsDevice_ = typename Base::WeightsDeviceView(
129 "rowWeightsDevice_", inmatrix->getLocalNumRows(),
130 this->nWeightsPerRow_);
131
132 this->numNzWeight_ = Kokkos::View<bool *, host_t>(
133 "numNzWeight_", this->nWeightsPerRow_);
134
135 for (int i = 0; i < this->nWeightsPerRow_; ++i) {
136 this->numNzWeight_(i) = false;
137 }
138 }
139}
140
142template <typename User, typename UserCoord>
143 template <typename Adapter>
145 const User &in, User *&out,
146 const PartitioningSolution<Adapter> &solution) const
147{
148 // Get an import list (rows to be received)
149 size_t numNewRows;
150 ArrayRCP<gno_t> importList;
151 try{
152 numNewRows = Zoltan2::getImportList<Adapter,
154 (solution, this, importList);
155 }
157
158 // Move the rows, creating a new matrix.
159 RCP<User> outPtr = this->doMigration(in, numNewRows,importList.getRawPtr());
160 out = const_cast<User *>(outPtr.get());
161 outPtr.release();
162}
163
165template <typename User, typename UserCoord>
166 template <typename Adapter>
168 const User &in, RCP<User> &out,
169 const PartitioningSolution<Adapter> &solution) const
170{
171 // Get an import list (rows to be received)
172 size_t numNewRows;
173 ArrayRCP<gno_t> importList;
174 try{
175 numNewRows = Zoltan2::getImportList<Adapter,
177 (solution, this, importList);
178 }
180
181 // Move the rows, creating a new matrix.
182 out = this->doMigration(in, numNewRows, importList.getRawPtr());
183}
184
185} //namespace Zoltan2
186
187#endif
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tmatrix_t
#define Z2_FORWARD_EXCEPTIONS
Forward an exception back through call stack.
Traits for application input objects.
Defines the MatrixAdapter interface.
Helper functions for Partitioning Problems.
This file defines the StridedData class.
Defines the TpetraRowMatrixAdapter class.
typename InputTraits< User >::node_t node_t
typename InputTraits< User >::lno_t lno_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.
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
void init(const RCP< const User > &inmatrix)
void applyPartitioningSolution(const User &in, User *&out, const PartitioningSolution< Adapter > &solution) const
RCP< const User > getUserMatrix() const
Access to user's matrix.
TpetraCrsMatrixAdapter(const RCP< const User > &inmatrix, int nWeightsPerRow=0)
Constructor.
Provides access for Zoltan2 to Tpetra::RowMatrix data.
Created by mbenlioglu on Aug 31, 2020.
size_t getImportList(const PartitioningSolution< SolutionAdapter > &solution, const DataAdapter *const data, ArrayRCP< typename DataAdapter::gno_t > &imports)
From a PartitioningSolution, get a list of IDs to be imported. Assumes part numbers in PartitioningSo...
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.