Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
TpetraCrsMatrixInput.cpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 //
46 // Basic testing of Zoltan2::TpetraCrsMatrixAdapter
47 
53 #include <Zoltan2_InputTraits.hpp>
54 #include <Zoltan2_TestHelpers.hpp>
57 
58 #include <Teuchos_Comm.hpp>
59 #include <Teuchos_CommHelpers.hpp>
60 #include <Teuchos_DefaultComm.hpp>
61 #include <Teuchos_RCP.hpp>
62 #include <cstdlib>
63 #include <stdexcept>
64 
65 using Teuchos::Comm;
66 using Teuchos::Comm;
67 using Teuchos::RCP;
68 using Teuchos::rcp;
69 using Teuchos::rcp_const_cast;
70 using Teuchos::rcp_dynamic_cast;
71 
72 using ztcrsmatrix_t = Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t>;
74 using device_t = typename node_t::device_type;
76 using execspace_t =
77  typename crsAdapter_t::ConstWeightsHostView1D::execution_space;
78 
80 
81 template <typename adapter_t, typename matrix_t>
82 void TestMatrixIds(adapter_t &ia, matrix_t &matrix) {
83 
84  using idsHost_t = typename adapter_t::ConstIdsHostView;
85  using offsetsHost_t = typename adapter_t::ConstOffsetsHostView;
86  using localInds_t =
87  typename adapter_t::user_t::nonconst_local_inds_host_view_type;
88  using localVals_t =
89  typename adapter_t::user_t::nonconst_values_host_view_type;
90 
91 
92  const auto nrows = matrix.getLocalNumRows();
93  const auto ncols = matrix.getLocalNumEntries();
94  const auto maxNumEntries = matrix.getLocalMaxNumRowEntries();
95 
96  typename adapter_t::Base::ConstIdsHostView colIdsHost_("colIdsHost_", ncols);
97  typename adapter_t::Base::ConstOffsetsHostView offsHost_("offsHost_",
98  nrows + 1);
99 
100  localInds_t localColInds("localColInds", maxNumEntries);
101  localVals_t localVals("localVals", maxNumEntries);
102 
103  for (size_t r = 0; r < nrows; r++) {
104  size_t numEntries = 0;
105  matrix.getLocalRowCopy(r, localColInds, localVals, numEntries);;
106 
107  offsHost_(r + 1) = offsHost_(r) + numEntries;
108  for (size_t e = offsHost_(r), i = 0; e < offsHost_(r + 1); e++) {
109  colIdsHost_(e) = matrix.getColMap()->getGlobalElement(localColInds(i++));
110  }
111  }
112 
113  idsHost_t rowIdsHost;
114  ia.getRowIDsHostView(rowIdsHost);
115 
116  const auto matrixIDS = matrix.getRowMap()->getLocalElementList();
117 
118  Z2_TEST_COMPARE_ARRAYS(matrixIDS, rowIdsHost);
119 
120  idsHost_t colIdsHost;
121  offsetsHost_t offsetsHost;
122  ia.getCRSHostView(offsetsHost, colIdsHost);
123 
124  Z2_TEST_COMPARE_ARRAYS(colIdsHost_, colIdsHost);
125  Z2_TEST_COMPARE_ARRAYS(offsHost_, offsetsHost);
126 }
127 
128 template <typename adapter_t, typename matrix_t>
129 void verifyInputAdapter(adapter_t &ia, matrix_t &matrix) {
130  using idsDevice_t = typename adapter_t::ConstIdsDeviceView;
131  using idsHost_t = typename adapter_t::ConstIdsHostView;
132  using weightsDevice_t = typename adapter_t::WeightsDeviceView1D;
133  using weightsHost_t = typename adapter_t::WeightsHostView1D;
134 
135  const auto nrows = ia.getLocalNumIDs();
136 
137  Z2_TEST_EQUALITY(ia.getLocalNumRows(), matrix.getLocalNumRows());
138  Z2_TEST_EQUALITY(ia.getLocalNumColumns(), matrix.getLocalNumCols());
139  Z2_TEST_EQUALITY(ia.getLocalNumEntries(), matrix.getLocalNumEntries());
140 
144 
145  idsDevice_t rowIdsDevice;
146  ia.getRowIDsDeviceView(rowIdsDevice);
147  idsHost_t rowIdsHost;
148  ia.getRowIDsHostView(rowIdsHost);
149 
150  Z2_TEST_DEVICE_HOST_VIEWS(rowIdsDevice, rowIdsHost);
151 
155  Z2_TEST_THROW(ia.setRowWeightsDevice(
156  typename adapter_t::WeightsDeviceView1D{}, 50),
157  std::runtime_error);
158 
159  weightsDevice_t wgts0("wgts0", nrows);
160  Kokkos::parallel_for(
161  nrows, KOKKOS_LAMBDA(const int idx) { wgts0(idx) = idx * 2; });
162 
163  Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts0, 0));
164 
165  // Don't reuse the same View, since we don't copy the values,
166  // we just assign the View (increase use count)
167  weightsDevice_t wgts1("wgts1", nrows);
168  Kokkos::parallel_for(
169  nrows, KOKKOS_LAMBDA(const int idx) { wgts1(idx) = idx * 3; });
170 
171  Z2_TEST_NOTHROW(ia.setRowWeightsDevice(wgts1, 1));
172 
176  {
177  weightsDevice_t weightsDevice;
178  Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 0));
179 
180  weightsHost_t weightsHost;
181  Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 0));
182 
183  Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
184 
185  Z2_TEST_DEVICE_HOST_VIEWS(wgts0, weightsHost);
186  }
187  {
188  weightsDevice_t weightsDevice;
189  Z2_TEST_NOTHROW(ia.getRowWeightsDeviceView(weightsDevice, 1));
190 
191  weightsHost_t weightsHost;
192  Z2_TEST_NOTHROW(ia.getRowWeightsHostView(weightsHost, 1));
193 
194  Z2_TEST_DEVICE_HOST_VIEWS(weightsDevice, weightsHost);
195 
196  Z2_TEST_DEVICE_HOST_VIEWS(wgts1, weightsHost);
197  }
198  {
199  weightsDevice_t wgtsDevice;
200  Z2_TEST_THROW(ia.getRowWeightsDeviceView(wgtsDevice, 2),
201  std::runtime_error);
202 
203  weightsHost_t wgtsHost;
204  Z2_TEST_THROW(ia.getRowWeightsHostView(wgtsHost, 2), std::runtime_error);
205  }
206 
207  TestMatrixIds(ia, matrix);
208 }
209 
211 
212 int main(int narg, char *arg[]) {
214  using crsPart_t = crsAdapter_t::part_t;
215 
216  Tpetra::ScopeGuard tscope(&narg, &arg);
217  const auto comm = Tpetra::getDefaultComm();
218 
219  try {
220  Teuchos::ParameterList params;
221  params.set("input file", "simple");
222  params.set("file type", "Chaco");
223 
224  auto uinput = rcp(new UserInputForTests(params, comm));
225 
226  // Input crs matrix
227  const auto crsMatrix = uinput->getUITpetraCrsMatrix();
228 
229  const auto nrows = crsMatrix->getLocalNumRows();
230 
231  // To test migration in the input adapter we need a Solution object.
232  const auto env = rcp(new Zoltan2::Environment(comm));
233 
234  const int nWeights = 2;
235 
237  // User object is Tpetra::CrsMatrix
239 
240  PrintFromRoot("Input adapter for Tpetra::CrsMatrix");
241 
242  // Graph Adapters use crsGraph, original TpetraInput uses trM (=CrsMatrix)
243  auto tpetraCrsMatrixInput = rcp(new crsAdapter_t(crsMatrix, nWeights));
244 
245  verifyInputAdapter(*tpetraCrsMatrixInput, *crsMatrix);
246 
247  crsPart_t *p = new crsPart_t[nrows];
248  memset(p, 0, sizeof(crsPart_t) * nrows);
249  ArrayRCP<crsPart_t> solnParts(p, 0, nrows, true);
250 
251  crsSoln_t solution(env, comm, nWeights);
252  solution.setParts(solnParts);
253 
254  ztcrsmatrix_t *mMigrate = NULL;
255  tpetraCrsMatrixInput->applyPartitioningSolution(*crsMatrix, mMigrate,
256  solution);
257  const auto newM = rcp(mMigrate);
258  auto cnewM = rcp_const_cast<const ztcrsmatrix_t>(newM);
259  auto newInput = rcp(new crsAdapter_t(cnewM, nWeights));
260 
261  PrintFromRoot("Input adapter for Tpetra::CrsMatrix migrated to proc 0");
262 
263  verifyInputAdapter(*newInput, *newM);
264 
265  } catch (std::exception &e) {
266  std::cout << e.what() << std::endl;
267  return EXIT_FAILURE;
268  }
269 
270  PrintFromRoot("PASS");
271 
272 }
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
void PrintFromRoot(const std::string &message)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t, znode_t > ztcrsmatrix_t
typename crsAdapter_t::ConstWeightsHostView1D::execution_space execspace_t
typename node_t::device_type device_t
Provides access for Zoltan2 to Tpetra::CrsMatrix data.
typename Zoltan2::InputTraits< ztcrsmatrix_t >::node_t node_t
int main(int narg, char **arg)
Definition: coloring1.cpp:199
common code used by tests
typename InputTraits< User >::part_t part_t
#define Z2_TEST_COMPARE_ARRAYS(val1, val2)
A PartitioningSolution is a solution to a partitioning problem.
#define Z2_TEST_DEVICE_HOST_VIEWS(deviceView, hostView)
#define Z2_TEST_EQUALITY(val1, val2)
Traits for application input objects.
Defines the TpetraRowMatrixAdapter class.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
Zoltan2::TpetraCrsMatrixAdapter< ztcrsmatrix_t > crsAdapter_t
default_node_t node_t
The Kokkos node type. This is only meaningful for users of Tpetra objects.
#define Z2_TEST_THROW(code, ExceptType)
void TestMatrixIds(adapter_t &ia, matrix_t &matrix)
#define Z2_TEST_NOTHROW(code)