Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
MatrixAdapter.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 // Testing of GraphAdapter built from Xpetra matrix input adapters.
12 //
13 
18 #include "Kokkos_StaticCrsGraph.hpp"
19 #include "Kokkos_UnorderedMap.hpp"
21 #include <Zoltan2_InputTraits.hpp>
22 #include <Zoltan2_TestHelpers.hpp>
25 
26 #include <bitset>
27 #include <iostream>
28 #include <string>
29 
30 #include <Teuchos_ArrayView.hpp>
31 #include <Teuchos_Comm.hpp>
32 #include <Teuchos_DefaultComm.hpp>
33 #include <Teuchos_TestingHelpers.hpp>
34 
35 using Teuchos::ArrayView;
36 using Teuchos::Comm;
37 using Teuchos::RCP;
38 using Teuchos::rcp;
39 using Teuchos::rcp_const_cast;
40 using Teuchos::rcp_dynamic_cast;
41 
43 
44 using tcrsMatrix_t = Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>;
45 using trowMatrix_t = Tpetra::RowMatrix<zscalar_t, zlno_t, zgno_t, znode_t>;
46 using tmap_t = Tpetra::Map<zlno_t, zgno_t, znode_t>;
47 
52 
54 int main(int narg, char *arg[]) {
55  Tpetra::ScopeGuard tscope(&narg, &arg);
56  Teuchos::RCP<const Teuchos::Comm<int>> comm = Tpetra::getDefaultComm();
57 
58  int rank = comm->getRank();
59 
60  int nVtxWeights = 2;
61  int nnzWgtIdx = -1;
62  std::string fname("simple");
63 
64  if (rank == 0)
65  std::cout << "TESTING base case (global)" << std::endl;
66 
67  // Input generator
68  UserInputForTests *uinput;
69 
70  uinput = new UserInputForTests(testDataFilePath, fname, comm, true);
71 
72  RCP<tcrsMatrix_t> M = uinput->getUITpetraCrsMatrix();
73  RCP<trowMatrix_t> trM = rcp_dynamic_cast<trowMatrix_t>(M);
74  RCP<const trowMatrix_t> ctrM = rcp_const_cast<const trowMatrix_t>(trM);
75  zlno_t nLocalRows = M->getLocalNumRows();
76  std::cout << "nLocalRows: " << nLocalRows << std::endl;
77 
78  // Weights:
79  zscalar_t **rowWeights = nullptr;
80  // create as many 1-D weights views as nVtxWeights
82  nLocalRows);
84  nLocalRows);
85 
86  auto wgts0Host = Kokkos::create_mirror_view(wgts0);
87  auto wgts1Host = Kokkos::create_mirror_view(wgts1);
88 
89  for (zlno_t i = 0; i < nLocalRows; i++) {
90  wgts0Host(i) = i;
91  wgts1Host(i) = 200000 + i;
92  }
93 
94  if (nVtxWeights > 0) {
95  rowWeights = new zscalar_t *[nVtxWeights];
96  for (int i = 0; i < nVtxWeights; i++) {
97  if (nnzWgtIdx == i) {
98  rowWeights[i] = nullptr;
99  } else {
100  rowWeights[i] = new zscalar_t[nLocalRows];
101  for (zlno_t j = 0; j < nLocalRows; j++) {
102  rowWeights[i][j] = 200000 * i + j;
103  }
104  }
105  }
106  }
107 
108  Kokkos::deep_copy(wgts0, wgts0Host);
109  Kokkos::deep_copy(wgts1, wgts1Host);
110 
111  tRowMAdapter_t tmi(ctrM, nVtxWeights);
112  for (int i = 0; i < nVtxWeights; i++) {
113  tmi.setWeights(rowWeights[i], 1, i);
114  }
115  tmi.setRowWeightsDevice(wgts0, 0);
116  tmi.setRowWeightsDevice(wgts1, 1);
117 
118  simpleVAdapter_t *via = nullptr;
119 
120  // Set up some fake input
121  zscalar_t **coords = nullptr;
122  int coordDim = 3;
123 
124  if (coordDim > 0) {
125  coords = new zscalar_t *[coordDim];
126  for (int i = 0; i < coordDim; i++) {
127  coords[i] = new zscalar_t[nLocalRows];
128  for (zlno_t j = 0; j < nLocalRows; j++) {
129  coords[i][j] = 100000 * i + j;
130  }
131  }
132  }
133 
134  zgno_t *gids = nullptr;
135  if (coordDim > 0) {
136  gids = new zgno_t[nLocalRows];
137  for (zlno_t i = 0; i < nLocalRows; i++)
138  gids[i] = M->getRowMap()->getGlobalElement(i);
139  via = new simpleVAdapter_t(nLocalRows, gids, coords[0],
140  (coordDim > 1 ? coords[1] : nullptr),
141  (coordDim > 2 ? coords[2] : nullptr), 1, 1, 1);
142  tmi.setCoordinateInput(via);
143  }
144 
145  // TEST of getIDsView, getIDsKokkosView, getRowIDsView, getRowIDsHostView and
146  // getRowIDsDeviceView
147 
148  const zgno_t *ids;
149  const zgno_t *rowIds;
153 
154  tmi.getIDsView(ids);
155  tmi.getRowIDsView(rowIds);
156  tmi.getIDsKokkosView(kIds);
157  auto kIdsHost = Kokkos::create_mirror_view(kIds);
158  Kokkos::deep_copy(kIdsHost, kIds);
159 
160  tmi.getRowIDsHostView(kHostIds);
161  tmi.getRowIDsDeviceView(kDeviceIds);
162 
163  auto kDeviceIdsHost = Kokkos::create_mirror_view(kDeviceIds);
164  Kokkos::deep_copy(kDeviceIdsHost, kDeviceIds);
165 
166  bool success = true;
167  for (size_t i = 0; i < tmi.getLocalNumIDs(); ++i) {
168  TEUCHOS_TEST_EQUALITY(ids[i], kIdsHost(i), std::cout, success);
169  TEUCHOS_TEST_EQUALITY(ids[i], rowIds[i], std::cout, success);
170  TEUCHOS_TEST_EQUALITY(kIdsHost(i), kHostIds(i), std::cout, success);
171  TEUCHOS_TEST_EQUALITY(kIdsHost(i), kDeviceIdsHost(i), std::cout, success);
172  }
173  TEST_FAIL_AND_EXIT(*comm, success, "ids != vertexIds != kIds", 1)
174 
175  // TEST of getCRSView, getCRSHostView and getCRSDeviceView
176  // const zoffset_t *offsets;
177  ArrayRCP<const zoffset_t> offsets;
178  ArrayRCP<const zgno_t> colIds;
183  tmi.getCRSView(offsets, colIds);
184  tmi.getCRSHostView(kHostOffsets, kHostColIds);
185  tmi.getCRSDeviceView(kDeviceOffsets, kDeviceColIds);
186 
187  auto kDeviceColIdsHost = Kokkos::create_mirror_view(kDeviceColIds);
188  Kokkos::deep_copy(kDeviceColIdsHost, kDeviceColIds);
189 
190  auto kDeviceOffsetsHost = Kokkos::create_mirror_view(kDeviceOffsets);
191  Kokkos::deep_copy(kDeviceOffsetsHost, kDeviceOffsets);
192 
193  for (int i = 0; success && i < colIds.size(); i++) {
194  TEUCHOS_TEST_EQUALITY(colIds[i], kHostColIds(i), std::cout, success);
195  TEUCHOS_TEST_EQUALITY(colIds[i], kDeviceColIdsHost(i), std::cout, success);
196  }
197  TEST_FAIL_AND_EXIT(*comm, success, "colIds != kHostColIds != kDeviceColIds",
198  1)
199 
200  for (int i = 0; success && i < offsets.size(); i++) {
201  TEUCHOS_TEST_EQUALITY(offsets[i], kHostOffsets(i), std::cout, success);
202  TEUCHOS_TEST_EQUALITY(offsets[i], kDeviceOffsetsHost(i), std::cout,
203  success);
204  }
205  TEST_FAIL_AND_EXIT(*comm, success,
206  "offsets != kHostOffsets != kDeviceOffsets", 1)
207 
208  ArrayRCP<const zscalar_t> values;
211 
212  tmi.getCRSView(offsets, colIds, values);
213  tmi.getCRSHostView(kHostOffsets, kHostColIds, kHostValues);
214  tmi.getCRSDeviceView(kDeviceOffsets, kDeviceColIds, kDeviceValues);
215  auto kDeviceValuesHost = Kokkos::create_mirror_view(kDeviceValues);
216  Kokkos::deep_copy(kDeviceValuesHost, kDeviceValues);
217 
218  for (int i = 0; success && i < colIds.size(); i++) {
219  TEUCHOS_TEST_EQUALITY(colIds[i], kHostColIds(i), std::cout, success);
220  TEUCHOS_TEST_EQUALITY(colIds[i], kDeviceColIdsHost(i), std::cout, success);
221  }
222  TEST_FAIL_AND_EXIT(*comm, success, "colIds != kHostColIds != kDeviceColIds",
223  1)
224 
225  for (int i = 0; success && i < offsets.size(); i++) {
226  TEUCHOS_TEST_EQUALITY(offsets[i], kHostOffsets(i), std::cout, success);
227  TEUCHOS_TEST_EQUALITY(offsets[i], kDeviceOffsetsHost(i), std::cout,
228  success);
229  }
230  TEST_FAIL_AND_EXIT(*comm, success,
231  "offsets != kHostOffsets != kDeviceOffsets", 1)
232 
233  for (int i = 0; success && i < values.size(); i++) {
234  TEUCHOS_TEST_EQUALITY(values[i], kHostValues(i), std::cout, success);
235  TEUCHOS_TEST_EQUALITY(values[i], kDeviceValuesHost(i), std::cout, success);
236  }
237  TEST_FAIL_AND_EXIT(*comm, success, "values != kHostValues != kDeviceValues",
238  1)
239 
240  // TEST of getRowWeightsView, getRowWeightsHost0View and
241  // getRowWeightsDeviceView
242 
247 
248  tmi.getRowWeightsHostView(weightsHost0, 0);
249  tmi.getRowWeightsDeviceView(weightsDevice0, 0);
250 
251  tmi.getRowWeightsHostView(weightsHost1, 1);
252  tmi.getRowWeightsDeviceView(weightsDevice1, 1);
253 
254  auto hostWeightsDevice0 = Kokkos::create_mirror_view(weightsDevice0);
255  Kokkos::deep_copy(hostWeightsDevice0, weightsDevice0);
256 
257  auto hostWeightsDevice1 = Kokkos::create_mirror_view(weightsDevice1);
258  Kokkos::deep_copy(hostWeightsDevice1, weightsDevice1);
259 
260  for (int w = 0; success && w < nVtxWeights; w++) {
261  const zscalar_t *wgts;
262 
263  Kokkos::View<zscalar_t *, typename znode_t::device_type> wkgts;
264  int stride;
265  tmi.getRowWeightsView(wgts, stride, w);
266 
267  if (w == 0) {
268 
269  for (zlno_t i = 0; success && i < nLocalRows; ++i) {
270  TEUCHOS_TEST_EQUALITY(wgts[stride * i], weightsHost0(i), std::cout,
271  success);
272  TEUCHOS_TEST_EQUALITY(wgts[stride * i], hostWeightsDevice0(i), std::cout,
273  success);
274  }
275  } else {
276  for (zlno_t i = 0; success && i < nLocalRows; ++i) {
277  TEUCHOS_TEST_EQUALITY(wgts[stride * i], weightsHost1(i), std::cout,
278  success);
279  TEUCHOS_TEST_EQUALITY(wgts[stride * i], hostWeightsDevice1(i), std::cout,
280  success);
281  }
282  }
283  }
284  TEST_FAIL_AND_EXIT(*comm, success, "wgts != vwgts != wkgts", 1)
285 
286  // TEST of useNumNonzerosAsRowWeight, useNumNonzerosAsRowWeightHost and
287  // useNumNonzerosAsRowWeightDevice
288  // for (int w = 0; success && w < nVtxWeights; w++) {
289  // TEUCHOS_TEST_EQUALITY(tmi.useNumNonzerosAsRowWeight(w),
290  // tmi.useNumNonzerosAsRowWeightHost(w), std::cout,
291  // success);
292  // TEUCHOS_TEST_EQUALITY(tmi.useNumNonzerosAsRowWeight(w),
293  // tmi.useNumNonzerosAsRowWeightDevice(w), std::cout,
294  // success);
295  // }
296  // TEST_FAIL_AND_EXIT(
297  // *comm, success,
298  // "useNumNonzerosAsRowWeight != useNumNonzerosAsRowWeightHost != "
299  // "useNumNonzerosAsRowWeightDevice",
300  // 1)
301 
302  // clean
303  if (nVtxWeights > 0) {
304  for (int i = 0; i < nVtxWeights; i++) {
305  if (rowWeights[i])
306  delete[] rowWeights[i];
307  }
308  delete[] rowWeights;
309  }
310 
311  if (coordDim > 0) {
312  delete via;
313  delete[] gids;
314  for (int i = 0; i < coordDim; i++) {
315  if (coords[i])
316  delete[] coords[i];
317  }
318  delete[] coords;
319  }
320 
321  delete uinput;
322 
323  if (rank == 0)
324  std::cout << "PASS" << std::endl;
325 
326  return 0;
327 }
typename baseMAdapter_t::offset_t zoffset_t
typename Zoltan2::TpetraRowMatrixAdapter< trowMatrix_t, simpleUser_t > tRowMAdapter_t
MatrixAdapter defines the adapter interface for matrices.
RCP< tcrsMatrix_t > getUITpetraCrsMatrix()
Kokkos::View< scalar_t *, device_t > WeightsDeviceView1D
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int main(int narg, char **arg)
Definition: coloring1.cpp:164
common code used by tests
Kokkos::View< const gno_t *, device_t > ConstIdsDeviceView
typename ConstScalarsDeviceView::HostMirror ConstScalarsHostView
Provides access for Zoltan2 to Tpetra::RowMatrix data.
list fname
Begin.
Definition: validXML.py:19
Defines the XpetraCrsMatrixAdapter class.
Kokkos::View< const scalar_t *, device_t > ConstScalarsDeviceView
Zoltan2::BasicVectorAdapter< simpleUser_t > simpleVAdapter_t
typename ConstOffsetsDeviceView::HostMirror ConstOffsetsHostView
typename ConstIdsDeviceView::HostMirror ConstIdsHostView
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Traits for application input objects.
Defines the TpetraRowMatrixAdapter class.
Tpetra::Map::local_ordinal_type zlno_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > tcrsMatrix_t
Tpetra::Map< zlno_t, zgno_t, znode_t > tmap_t
typename InputTraits< User >::offset_t offset_t
typename WeightsDeviceView1D::HostMirror WeightsHostView1D
float zscalar_t
Defines the BasicVectorAdapter class.
Kokkos::View< const offset_t *, device_t > ConstOffsetsDeviceView
Tpetra::Map::global_ordinal_type zgno_t
Tpetra::RowMatrix< zscalar_t, zlno_t, zgno_t, znode_t > trowMatrix_t
std::string testDataFilePath(".")