Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
BasicVectorAdapter.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 // Test for Zoltan2::BasicVectorAdapter
47 
49 #include <Zoltan2_TestHelpers.hpp>
50 
51 #include <Teuchos_DefaultComm.hpp>
52 #include <Teuchos_RCP.hpp>
53 #include <Teuchos_CommHelpers.hpp>
54 #include <Teuchos_TestingHelpers.hpp>
55 
56 #include <iostream>
57 
59 
60 void testBasisVector(Zoltan2::BasicVectorAdapter<userTypes_t> *ia, int *valueStrides, int wdim, int mvdim) {
61 
62  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
63  // Ids
64  const zgno_t *ids;
67 
68  ia->getIDsView(ids);
69  ia->getIDsHostView(kHostIds);
70  ia->getIDsDeviceView(kDeviceIds);
71 
72  auto kDeviceIdsHost = Kokkos::create_mirror_view(kDeviceIds);
73  Kokkos::deep_copy(kDeviceIdsHost, kDeviceIds);
74 
75  bool success = true;
76  for (size_t i = 0; i < ia->getLocalNumIDs(); ++i) {
77  TEUCHOS_TEST_EQUALITY(ids[i], kHostIds(i), std::cout, success);
78  TEUCHOS_TEST_EQUALITY(kHostIds(i), kDeviceIdsHost(i), std::cout, success);
79  }
80  TEST_FAIL_AND_EXIT(*comm, success, "ids != hostIds != deviceIds", 1)
81 
82  // Weights
85  ia->getWeightsHostView(kHostWgts);
86  ia->getWeightsDeviceView(kDeviceWgts);
87  auto kDeviceWgtsHost = Kokkos::create_mirror_view(kDeviceWgts);
88  Kokkos::deep_copy(kDeviceWgtsHost, kDeviceWgts);
89 
90  for (int w = 0; success && w < wdim; w++) {
91  const zscalar_t *wgts;
92 
93  Kokkos::View<zscalar_t **, typename znode_t::device_type> wkgts;
94  int stride;
95  ia->getWeightsView(wgts, stride, w);
96  for (zlno_t i = 0; success && i < ia->getNumWeightsPerID(); ++i) {
97  TEUCHOS_TEST_EQUALITY(wgts[stride * i], kHostWgts(i, w), std::cout,
98  success);
99  TEUCHOS_TEST_EQUALITY(wgts[stride * i], kDeviceWgtsHost(i, w), std::cout,
100  success);
101  }
102  }
103  TEST_FAIL_AND_EXIT(*comm, success, "wgts != vwgts != wkgts", 1)
104 
105  // Coords
108  ia->getCoordinatesHostView(kHostCoords);
109  ia->getCoordinatesDeviceView(kDeviceCoords);
110  auto kHostCoordsMV = Kokkos::create_mirror_view(kHostCoords);
111  Kokkos::deep_copy(kHostCoordsMV, kHostCoords);
112 
113  auto kDeviceCoordsMV = Kokkos::create_mirror_view(kDeviceCoords);
114  Kokkos::deep_copy(kDeviceCoordsMV, kDeviceCoords);
115  for (int v=0; v < mvdim; v++){
116  const zscalar_t *coords;
117  int stride;
118 
119  ia->getCoordinatesView(coords, stride, v);
120 
121  success = true;
122  for (size_t i = 0; i < ia->getLocalNumIDs(); ++i) {
123  TEUCHOS_TEST_EQUALITY(coords[i*stride], kHostCoordsMV(i, v), std::cout, success);
124  TEUCHOS_TEST_EQUALITY(coords[i*stride], kDeviceCoordsMV(i, v), std::cout, success);
125  }
126  TEST_FAIL_AND_EXIT(*comm, success, "ids != kHostCoordsMV != kDeviceCoordsMV", 1)
127  }
128 }
129 
130 int main(int narg, char *arg[])
131 {
132  Tpetra::ScopeGuard tscope(&narg, &arg);
133  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
134 
135  int rank = comm->getRank();
136  int nprocs = comm->getSize();
137  int fail = 0;
138 
139  // Create a single vector and a strided multi-vector with
140  // strided multi-weights.
141 
142  zlno_t numLocalIds = 10;
143  zgno_t *myIds = new zgno_t[numLocalIds];
144  zgno_t myFirstId = rank * numLocalIds;
145 
146  int wdim = 2;
147  zscalar_t *weights = new zscalar_t [numLocalIds*wdim];
148  int *weightStrides = new int [wdim];
149  const zscalar_t **weightPtrs = new const zscalar_t * [wdim];
150 
151  int mvdim = 3;
152  zscalar_t *v_values= new zscalar_t [numLocalIds];
153  zscalar_t *mv_values= new zscalar_t [mvdim * numLocalIds];
154  int *valueStrides = new int [mvdim];
155  const zscalar_t **valuePtrs = new const zscalar_t * [mvdim];
156 
157  for (zlno_t i=0; i < numLocalIds; i++){
158  myIds[i] = myFirstId+i;
159 
160  for (int w=0; w < wdim; w++)
161  weights[w*numLocalIds + i] = w + 1 + nprocs - rank;
162 
163  v_values[i] = numLocalIds-i;
164 
165  for (int v=0; v < mvdim; v++)
166  mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1);
167  }
168 
169  for (int w=0; w < wdim; w++){
170  weightStrides[w] = 1;
171  weightPtrs[w] = weights + numLocalIds*w;
172  }
173 
174  for (int v=0; v < mvdim; v++){
175  valueStrides[v] = mvdim;
176  valuePtrs[v] = mv_values + v;
177  }
178 
180 
181  {
182  // A Zoltan2::BasicVectorAdapter object with one vector and weights
183 
184  std::vector<const zscalar_t *> weightValues;
185  std::vector<int> strides;
186 
187  weightValues.push_back(weightPtrs[0]);
188  weightValues.push_back(weightPtrs[1]);
189  strides.push_back(1);
190  strides.push_back(1);
191 
192  try{
193  ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
194  v_values, 1, true, weightPtrs[0], 1);
195  }
196  catch (std::exception &e){
197  fail = 1;
198  }
199 
200  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
201 
202  testBasisVector(ia, valueStrides, 1, 1);
203 
204  delete ia;
205 
206 
207  }
208 
209  {
210  // A Zoltan2::BasicVectorAdapter object with a multivector with weights
211 
212  std::vector<const zscalar_t *> weightValues, values;
213  std::vector<int> wstrides, vstrides;
214 
215  for (int dim=0; dim < wdim; dim++){
216  weightValues.push_back(weightPtrs[dim]);
217  wstrides.push_back(1);
218  }
219 
220  for (int dim=0; dim < mvdim; dim++){
221  values.push_back(valuePtrs[dim]);
222  vstrides.push_back(mvdim);
223  }
224 
225  try{
227  numLocalIds, myIds, values, vstrides, weightValues, wstrides);
228  }
229  catch (std::exception &e){
230  fail = 1;
231  }
232 
233  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
234 
235  testBasisVector(ia, valueStrides, wdim, mvdim);
236 
237  delete ia;
238 
239  }
240 
241  if (rank == 0)
242  std::cout << "PASS" << std::endl;
243 
244  delete [] myIds;
245  delete [] weights;
246  delete [] weightStrides;
247  delete [] weightPtrs;
248  delete [] v_values;
249  delete [] mv_values;
250  delete [] valueStrides;
251  delete [] valuePtrs;
252 
253  return fail;
254 }
255 
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getIDsHostView(typename Base::ConstIdsHostView &ids) const override
void getCoordinatesView(const scalar_t *&elements, int &stride, int idx=0) const override
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int main(int narg, char **arg)
Definition: coloring1.cpp:199
void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
common code used by tests
void getIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override
Kokkos::View< const gno_t *, device_t > ConstIdsDeviceView
typename WeightsDeviceView::HostMirror WeightsHostView
typename CoordsDeviceView::HostMirror CoordsHostView
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
typename ConstIdsDeviceView::HostMirror ConstIdsHostView
void getCoordinatesHostView(typename AdapterWithCoords< User >::CoordsHostView &elements) const override
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Tpetra::Map::local_ordinal_type zlno_t
static const std::string fail
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process&#39; identifiers.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
void testBasisVector(Zoltan2::BasicVectorAdapter< userTypes_t > *ia, int *valueStrides, int wdim, int mvdim)
Kokkos::View< scalar_t **, device_t > WeightsDeviceView
float zscalar_t
Defines the BasicVectorAdapter class.
Tpetra::Map::global_ordinal_type zgno_t
void getCoordinatesDeviceView(typename AdapterWithCoords< User >::CoordsDeviceView &elements) const override