Zoltan2
Loading...
Searching...
No Matches
TpetraRowMatrixInput.cpp
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
10//
11// Basic testing of Zoltan2::TpetraRowMatrixAdapter
12
22
23#include <Teuchos_Comm.hpp>
24#include <Teuchos_CommHelpers.hpp>
25#include <Teuchos_DefaultComm.hpp>
26#include <Teuchos_RCP.hpp>
27#include <cstdlib>
28#include <stdexcept>
29
30using Teuchos::Comm;
31using Teuchos::Comm;
32using Teuchos::RCP;
33using Teuchos::rcp;
34using Teuchos::rcp_const_cast;
35using Teuchos::rcp_dynamic_cast;
36
37using ztcrsmatrix_t = Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>;
38using ztrowmatrix_t = Tpetra::RowMatrix<zscalar_t, zlno_t, zgno_t, znode_t>;
40using device_t = typename node_t::device_type;
44 typename rowAdapter_t::ConstWeightsHostView1D::execution_space;
45
47
48template<typename offset_t>
49void printMatrix(RCP<const Comm<int> > &comm, zlno_t nrows,
50 const zgno_t *rowIds, const offset_t *offsets, const zgno_t *colIds) {
51 int rank = comm->getRank();
52 int nprocs = comm->getSize();
53 comm->barrier();
54 for (int p=0; p < nprocs; p++){
55 if (p == rank){
56 std::cout << rank << ":" << std::endl;
57 for (zlno_t i=0; i < nrows; i++){
58 std::cout << " row " << rowIds[i] << ": ";
59 for (offset_t j=offsets[i]; j < offsets[i+1]; j++){
60 std::cout << colIds[j] << " ";
61 }
62 std::cout << std::endl;
63 }
64 std::cout.flush();
65 }
66 comm->barrier();
67 }
68 comm->barrier();
69}
70
72
73template <typename adapter_t, typename matrix_t>
74void TestMatrixIds(adapter_t &ia, matrix_t &matrix) {
75
76 using idsHost_t = typename adapter_t::ConstIdsHostView;
77 using offsetsHost_t = typename adapter_t::ConstOffsetsHostView;
78 using localInds_t =
79 typename adapter_t::user_t::nonconst_local_inds_host_view_type;
80 using localVals_t =
81 typename adapter_t::user_t::nonconst_values_host_view_type;
82
83
84 const auto nrows = matrix.getLocalNumRows();
85 const auto ncols = matrix.getLocalNumEntries();
86 const auto maxNumEntries = matrix.getLocalMaxNumRowEntries();
87
88 typename adapter_t::Base::ConstIdsHostView colIdsHost_("colIdsHost_", ncols);
89 typename adapter_t::Base::ConstOffsetsHostView offsHost_("offsHost_",
90 nrows + 1);
91
92 localInds_t localColInds("localColInds", maxNumEntries);
93 localVals_t localVals("localVals", maxNumEntries);
94
95 for (size_t r = 0; r < nrows; r++) {
96 size_t numEntries = 0;
97 matrix.getLocalRowCopy(r, localColInds, localVals, numEntries);;
98
99 offsHost_(r + 1) = offsHost_(r) + numEntries;
100 for (size_t e = offsHost_(r), i = 0; e < offsHost_(r + 1); e++) {
101 colIdsHost_(e) = matrix.getColMap()->getGlobalElement(localColInds(i++));
102 }
103 }
104
105 idsHost_t rowIdsHost;
106 ia.getRowIDsHostView(rowIdsHost);
107
108 const auto matrixIDS = matrix.getRowMap()->getLocalElementList();
109
110 Z2_TEST_COMPARE_ARRAYS(matrixIDS, rowIdsHost);
111
112 idsHost_t colIdsHost;
113 offsetsHost_t offsetsHost;
114 ia.getCRSHostView(offsetsHost, colIdsHost);
115
116 Z2_TEST_COMPARE_ARRAYS(colIdsHost_, colIdsHost);
117 Z2_TEST_COMPARE_ARRAYS(offsHost_, offsetsHost);
118}
119
120template <typename adapter_t, typename matrix_t>
121void verifyInputAdapter(adapter_t &ia, matrix_t &matrix) {
122 using idsDevice_t = typename adapter_t::ConstIdsDeviceView;
123 using idsHost_t = typename adapter_t::ConstIdsHostView;
124 using weightsDevice_t = typename adapter_t::WeightsDeviceView1D;
125 using weightsHost_t = typename adapter_t::WeightsHostView1D;
126
127 const auto nrows = ia.getLocalNumIDs();
128
129 Z2_TEST_EQUALITY(ia.getLocalNumRows(), matrix.getLocalNumRows());
130 Z2_TEST_EQUALITY(ia.getLocalNumColumns(), matrix.getLocalNumCols());
131 Z2_TEST_EQUALITY(ia.getLocalNumEntries(), matrix.getLocalNumEntries());
132
136
137 idsDevice_t rowIdsDevice;
138 ia.getRowIDsDeviceView(rowIdsDevice);
139 idsHost_t rowIdsHost;
140 ia.getRowIDsHostView(rowIdsHost);
141
142 Z2_TEST_DEVICE_HOST_VIEWS(rowIdsDevice, rowIdsHost);
143
147 Z2_TEST_THROW(ia.setRowWeightsDevice(
148 typename adapter_t::WeightsDeviceView1D{}, 50),
149 std::runtime_error);
150
151 weightsDevice_t wgts0("wgts0", nrows);
152 Kokkos::parallel_for(
153 nrows, KOKKOS_LAMBDA(const int idx) { wgts0(idx) = idx * 2; });
154
155 Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts0, 0));
156
157 // Don't reuse the same View, since we don't copy the values,
158 // we just assign the View (increase use count)
159 weightsDevice_t wgts1("wgts1", nrows);
160 Kokkos::parallel_for(
161 nrows, KOKKOS_LAMBDA(const int idx) { wgts1(idx) = idx * 3; });
162
163 Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts1, 1));
164
168 {
169 weightsDevice_t weightsDevice;
170 Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 0));
171
172 weightsHost_t weightsHost;
173 Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 0));
174
175 Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
176
177 Z2_TEST_DEVICE_HOST_VIEWS(wgts0, weightsHost);
178 }
179 {
180 weightsDevice_t weightsDevice;
181 Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 1));
182
183 weightsHost_t weightsHost;
184 Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 1));
185
186 Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
187
188 Z2_TEST_DEVICE_HOST_VIEWS(wgts1, weightsHost);
189 }
190 {
191 weightsDevice_t wgtsDevice;
192 Z2_TEST_THROW(ia.getRowWeightsDeviceView(wgtsDevice, 2),
193 std::runtime_error);
194
195 weightsHost_t wgtsHost;
196 Z2_TEST_THROW(ia.getRowWeightsHostView(wgtsHost, 2), std::runtime_error);
197 }
198
199 TestMatrixIds(ia, matrix);
200}
201
203
204int main(int narg, char *arg[]) {
206 using rowPart_t = rowAdapter_t::part_t;
207
208 Tpetra::ScopeGuard tscope(&narg, &arg);
209 const auto comm = Tpetra::getDefaultComm();
210
211 try {
212 Teuchos::ParameterList params;
213 params.set("input file", "simple");
214 params.set("file type", "Chaco");
215
216 auto uinput = rcp(new UserInputForTests(params, comm));
217
218 // Input crs matrix and row matrix cast from it.
219 const auto crsMatrix = uinput->getUITpetraCrsMatrix();
220 const auto rowMatrix = rcp_dynamic_cast<ztrowmatrix_t>(crsMatrix);
221
222 const auto nrows = rowMatrix->getLocalNumRows();
223
224 // To test migration in the input adapter we need a Solution object.
225 const auto env = rcp(new Zoltan2::Environment(comm));
226
227 const int nWeights = 2;
228
230 // User object is Tpetra::RowMatrix
232
233 PrintFromRoot("Input adapter for Tpetra::RowMatrix");
234
235 // Graph Adapters use crsGraph, original TpetraInput uses trM (=rowMatrix)
236 auto tpetraRowMatrixInput = rcp(new rowAdapter_t(rowMatrix, nWeights));
237
238 verifyInputAdapter(*tpetraRowMatrixInput, *rowMatrix);
239
240 rowPart_t *p = new rowPart_t[nrows];
241 memset(p, 0, sizeof(rowPart_t) * nrows);
242 ArrayRCP<rowPart_t> solnParts(p, 0, nrows, true);
243
244 rowSoln_t solution(env, comm, nWeights);
245 solution.setParts(solnParts);
246
247 ztrowmatrix_t *mMigrate = NULL;
248 tpetraRowMatrixInput->applyPartitioningSolution(*rowMatrix, mMigrate,
249 solution);
250 const auto newM = rcp(mMigrate);
251 auto cnewM = rcp_const_cast<const ztrowmatrix_t>(newM);
252 auto newInput = rcp(new rowAdapter_t(cnewM, nWeights));
253
254 PrintFromRoot("Input adapter for Tpetra::RowMatrix migrated to proc 0");
255
256 verifyInputAdapter(*newInput, *newM);
257
258 } catch (std::exception &e) {
259 std::cout << e.what() << std::endl;
260 return EXIT_FAILURE;
261 }
262
263 PrintFromRoot("PASS");
264
265}
typename node_t::device_type device_t
typename crsAdapter_t::ConstWeightsHostView1D::execution_space execspace_t
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztcrsmatrix_t
Tpetra::RowMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztrowmatrix_t
void TestMatrixIds(adapter_t &ia, matrix_t &matrix)
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
void printMatrix(RCP< const Comm< int > > &comm, zlno_t nrows, const zgno_t *rowIds, const offset_t *offsets, const zgno_t *colIds)
Zoltan2::TpetraRowMatrixAdapter< ztrowmatrix_t > rowAdapter_t
Traits for application input objects.
common code used by tests
void PrintFromRoot(const std::string &message)
Tpetra::Map ::local_ordinal_type zlno_t
#define Z2_TEST_DEVICE_HOST_VIEWS(deviceView, hostView)
#define Z2_TEST_THROW(code, ExceptType)
#define Z2_TEST_COMPARE_ARRAYS(val1, val2)
#define Z2_TEST_EQUALITY(val1, val2)
Tpetra::Map ::global_ordinal_type zgno_t
#define Z2_TEST_NOTHROW(code)
Defines the TpetraRowMatrixAdapter class.
int main()
typename InputTraits< User >::part_t part_t
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
Provides access for Zoltan2 to Tpetra::RowGraph data.
Provides access for Zoltan2 to Tpetra::RowMatrix data.
size_t getLocalNumRows() const
Returns the number of rows on this process.
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.