Zoltan2
Loading...
Searching...
No Matches
TpetraCrsMatrixInput.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::TpetraCrsMatrixAdapter
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>;
39using device_t = typename node_t::device_type;
42 typename crsAdapter_t::ConstWeightsHostView1D::execution_space;
43
45
46template <typename adapter_t, typename matrix_t>
47void TestMatrixIds(adapter_t &ia, matrix_t &matrix) {
48
49 using idsHost_t = typename adapter_t::ConstIdsHostView;
50 using offsetsHost_t = typename adapter_t::ConstOffsetsHostView;
51 using localInds_t =
52 typename adapter_t::user_t::nonconst_local_inds_host_view_type;
53 using localVals_t =
54 typename adapter_t::user_t::nonconst_values_host_view_type;
55
56
57 const auto nrows = matrix.getLocalNumRows();
58 const auto ncols = matrix.getLocalNumEntries();
59 const auto maxNumEntries = matrix.getLocalMaxNumRowEntries();
60
61 typename adapter_t::Base::ConstIdsHostView colIdsHost_("colIdsHost_", ncols);
62 typename adapter_t::Base::ConstOffsetsHostView offsHost_("offsHost_",
63 nrows + 1);
64
65 localInds_t localColInds("localColInds", maxNumEntries);
66 localVals_t localVals("localVals", maxNumEntries);
67
68 for (size_t r = 0; r < nrows; r++) {
69 size_t numEntries = 0;
70 matrix.getLocalRowCopy(r, localColInds, localVals, numEntries);;
71
72 offsHost_(r + 1) = offsHost_(r) + numEntries;
73 for (size_t e = offsHost_(r), i = 0; e < offsHost_(r + 1); e++) {
74 colIdsHost_(e) = matrix.getColMap()->getGlobalElement(localColInds(i++));
75 }
76 }
77
78 idsHost_t rowIdsHost;
79 ia.getRowIDsHostView(rowIdsHost);
80
81 const auto matrixIDS = matrix.getRowMap()->getLocalElementList();
82
83 Z2_TEST_COMPARE_ARRAYS(matrixIDS, rowIdsHost);
84
85 idsHost_t colIdsHost;
86 offsetsHost_t offsetsHost;
87 ia.getCRSHostView(offsetsHost, colIdsHost);
88
89 Z2_TEST_COMPARE_ARRAYS(colIdsHost_, colIdsHost);
90 Z2_TEST_COMPARE_ARRAYS(offsHost_, offsetsHost);
91}
92
93template <typename adapter_t, typename matrix_t>
94void verifyInputAdapter(adapter_t &ia, matrix_t &matrix) {
95 using idsDevice_t = typename adapter_t::ConstIdsDeviceView;
96 using idsHost_t = typename adapter_t::ConstIdsHostView;
97 using weightsDevice_t = typename adapter_t::WeightsDeviceView1D;
98 using weightsHost_t = typename adapter_t::WeightsHostView1D;
99
100 const auto nrows = ia.getLocalNumIDs();
101
102 Z2_TEST_EQUALITY(ia.getLocalNumRows(), matrix.getLocalNumRows());
103 Z2_TEST_EQUALITY(ia.getLocalNumColumns(), matrix.getLocalNumCols());
104 Z2_TEST_EQUALITY(ia.getLocalNumEntries(), matrix.getLocalNumEntries());
105
109
110 idsDevice_t rowIdsDevice;
111 ia.getRowIDsDeviceView(rowIdsDevice);
112 idsHost_t rowIdsHost;
113 ia.getRowIDsHostView(rowIdsHost);
114
115 Z2_TEST_DEVICE_HOST_VIEWS(rowIdsDevice, rowIdsHost);
116
120 Z2_TEST_THROW(ia.setRowWeightsDevice(
121 typename adapter_t::WeightsDeviceView1D{}, 50),
122 std::runtime_error);
123
124 weightsDevice_t wgts0("wgts0", nrows);
125 Kokkos::parallel_for(
126 nrows, KOKKOS_LAMBDA(const int idx) { wgts0(idx) = idx * 2; });
127
128 Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts0, 0));
129
130 // Don't reuse the same View, since we don't copy the values,
131 // we just assign the View (increase use count)
132 weightsDevice_t wgts1("wgts1", nrows);
133 Kokkos::parallel_for(
134 nrows, KOKKOS_LAMBDA(const int idx) { wgts1(idx) = idx * 3; });
135
136 Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts1, 1));
137
141 {
142 weightsDevice_t weightsDevice;
143 Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 0));
144
145 weightsHost_t weightsHost;
146 Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 0));
147
148 Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
149
150 Z2_TEST_DEVICE_HOST_VIEWS(wgts0, weightsHost);
151 }
152 {
153 weightsDevice_t weightsDevice;
154 Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 1));
155
156 weightsHost_t weightsHost;
157 Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 1));
158
159 Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
160
161 Z2_TEST_DEVICE_HOST_VIEWS(wgts1, weightsHost);
162 }
163 {
164 weightsDevice_t wgtsDevice;
165 Z2_TEST_THROW(ia.getRowWeightsDeviceView(wgtsDevice, 2),
166 std::runtime_error);
167
168 weightsHost_t wgtsHost;
169 Z2_TEST_THROW(ia.getRowWeightsHostView(wgtsHost, 2), std::runtime_error);
170 }
171
172 TestMatrixIds(ia, matrix);
173}
174
176
177int main(int narg, char *arg[]) {
179 using crsPart_t = crsAdapter_t::part_t;
180
181 Tpetra::ScopeGuard tscope(&narg, &arg);
182 const auto comm = Tpetra::getDefaultComm();
183
184 try {
185 Teuchos::ParameterList params;
186 params.set("input file", "simple");
187 params.set("file type", "Chaco");
188
189 auto uinput = rcp(new UserInputForTests(params, comm));
190
191 // Input crs matrix
192 const auto crsMatrix = uinput->getUITpetraCrsMatrix();
193
194 const auto nrows = crsMatrix->getLocalNumRows();
195
196 // To test migration in the input adapter we need a Solution object.
197 const auto env = rcp(new Zoltan2::Environment(comm));
198
199 const int nWeights = 2;
200
202 // User object is Tpetra::CrsMatrix
204
205 PrintFromRoot("Input adapter for Tpetra::CrsMatrix");
206
207 // Graph Adapters use crsGraph, original TpetraInput uses trM (=CrsMatrix)
208 auto tpetraCrsMatrixInput = rcp(new crsAdapter_t(crsMatrix, nWeights));
209
210 verifyInputAdapter(*tpetraCrsMatrixInput, *crsMatrix);
211
212 crsPart_t *p = new crsPart_t[nrows];
213 memset(p, 0, sizeof(crsPart_t) * nrows);
214 ArrayRCP<crsPart_t> solnParts(p, 0, nrows, true);
215
216 crsSoln_t solution(env, comm, nWeights);
217 solution.setParts(solnParts);
218
219 ztcrsmatrix_t *mMigrate = NULL;
220 tpetraCrsMatrixInput->applyPartitioningSolution(*crsMatrix, mMigrate,
221 solution);
222 const auto newM = rcp(mMigrate);
223 auto cnewM = rcp_const_cast<const ztcrsmatrix_t>(newM);
224 auto newInput = rcp(new crsAdapter_t(cnewM, nWeights));
225
226 PrintFromRoot("Input adapter for Tpetra::CrsMatrix migrated to proc 0");
227
228 verifyInputAdapter(*newInput, *newM);
229
230 } catch (std::exception &e) {
231 std::cout << e.what() << std::endl;
232 return EXIT_FAILURE;
233 }
234
235 PrintFromRoot("PASS");
236
237}
Zoltan2::TpetraCrsMatrixAdapter< ztcrsmatrix_t > crsAdapter_t
typename node_t::device_type device_t
typename crsAdapter_t::ConstWeightsHostView1D::execution_space execspace_t
void TestMatrixIds(adapter_t &ia, matrix_t &matrix)
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztcrsmatrix_t
Traits for application input objects.
common code used by tests
void PrintFromRoot(const std::string &message)
#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)
#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.
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.