Zoltan2
Loading...
Searching...
No Matches
BasicVectorAdapter.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// Test for Zoltan2::BasicVectorAdapter
12
15
16#include <Teuchos_DefaultComm.hpp>
17#include <Teuchos_RCP.hpp>
18#include <Teuchos_CommHelpers.hpp>
19#include <Teuchos_TestingHelpers.hpp>
20
21#include <iostream>
22
24
25void testBasisVector(Zoltan2::BasicVectorAdapter<userTypes_t> *ia, int *valueStrides, int wdim, int mvdim) {
26
27 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
28 // Ids
29 const zgno_t *ids;
32
33 ia->getIDsView(ids);
34 ia->getIDsHostView(kHostIds);
35 ia->getIDsDeviceView(kDeviceIds);
36
37 auto kDeviceIdsHost = Kokkos::create_mirror_view(kDeviceIds);
38 Kokkos::deep_copy(kDeviceIdsHost, kDeviceIds);
39
40 bool success = true;
41 for (size_t i = 0; i < ia->getLocalNumIDs(); ++i) {
42 TEUCHOS_TEST_EQUALITY(ids[i], kHostIds(i), std::cout, success);
43 TEUCHOS_TEST_EQUALITY(kHostIds(i), kDeviceIdsHost(i), std::cout, success);
44 }
45 TEST_FAIL_AND_EXIT(*comm, success, "ids != hostIds != deviceIds", 1)
46
47 // Weights
50 ia->getWeightsHostView(kHostWgts);
51 ia->getWeightsDeviceView(kDeviceWgts);
52 auto kDeviceWgtsHost = Kokkos::create_mirror_view(kDeviceWgts);
53 Kokkos::deep_copy(kDeviceWgtsHost, kDeviceWgts);
54
55 for (int w = 0; success && w < wdim; w++) {
56 const zscalar_t *wgts;
57
58 Kokkos::View<zscalar_t **, typename znode_t::device_type> wkgts;
59 int stride;
60 ia->getWeightsView(wgts, stride, w);
61 for (zlno_t i = 0; success && i < ia->getNumWeightsPerID(); ++i) {
62 TEUCHOS_TEST_EQUALITY(wgts[stride * i], kHostWgts(i, w), std::cout,
63 success);
64 TEUCHOS_TEST_EQUALITY(wgts[stride * i], kDeviceWgtsHost(i, w), std::cout,
65 success);
66 }
67 }
68 TEST_FAIL_AND_EXIT(*comm, success, "wgts != vwgts != wkgts", 1)
69
70 // Coords
73 ia->getCoordinatesHostView(kHostCoords);
74 ia->getCoordinatesDeviceView(kDeviceCoords);
75 auto kHostCoordsMV = Kokkos::create_mirror_view(kHostCoords);
76 Kokkos::deep_copy(kHostCoordsMV, kHostCoords);
77
78 auto kDeviceCoordsMV = Kokkos::create_mirror_view(kDeviceCoords);
79 Kokkos::deep_copy(kDeviceCoordsMV, kDeviceCoords);
80 for (int v=0; v < mvdim; v++){
81 const zscalar_t *coords;
82 int stride;
83
84 ia->getCoordinatesView(coords, stride, v);
85
86 success = true;
87 for (size_t i = 0; i < ia->getLocalNumIDs(); ++i) {
88 TEUCHOS_TEST_EQUALITY(coords[i*stride], kHostCoordsMV(i, v), std::cout, success);
89 TEUCHOS_TEST_EQUALITY(coords[i*stride], kDeviceCoordsMV(i, v), std::cout, success);
90 }
91 TEST_FAIL_AND_EXIT(*comm, success, "ids != kHostCoordsMV != kDeviceCoordsMV", 1)
92 }
93}
94
95int main(int narg, char *arg[])
96{
97 Tpetra::ScopeGuard tscope(&narg, &arg);
98 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
99
100 int rank = comm->getRank();
101 int nprocs = comm->getSize();
102 int fail = 0;
103
104 // Create a single vector and a strided multi-vector with
105 // strided multi-weights.
106
107 zlno_t numLocalIds = 10;
108 zgno_t *myIds = new zgno_t[numLocalIds];
109 zgno_t myFirstId = rank * numLocalIds;
110
111 int wdim = 2;
112 zscalar_t *weights = new zscalar_t [numLocalIds*wdim];
113 int *weightStrides = new int [wdim];
114 const zscalar_t **weightPtrs = new const zscalar_t * [wdim];
115
116 int mvdim = 3;
117 zscalar_t *v_values= new zscalar_t [numLocalIds];
118 zscalar_t *mv_values= new zscalar_t [mvdim * numLocalIds];
119 int *valueStrides = new int [mvdim];
120 const zscalar_t **valuePtrs = new const zscalar_t * [mvdim];
121
122 for (zlno_t i=0; i < numLocalIds; i++){
123 myIds[i] = myFirstId+i;
124
125 for (int w=0; w < wdim; w++)
126 weights[w*numLocalIds + i] = w + 1 + nprocs - rank;
127
128 v_values[i] = numLocalIds-i;
129
130 for (int v=0; v < mvdim; v++)
131 mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1);
132 }
133
134 for (int w=0; w < wdim; w++){
135 weightStrides[w] = 1;
136 weightPtrs[w] = weights + numLocalIds*w;
137 }
138
139 for (int v=0; v < mvdim; v++){
140 valueStrides[v] = mvdim;
141 valuePtrs[v] = mv_values + v;
142 }
143
145
146 {
147 // A Zoltan2::BasicVectorAdapter object with one vector and weights
148
149 std::vector<const zscalar_t *> weightValues;
150 std::vector<int> strides;
151
152 weightValues.push_back(weightPtrs[0]);
153 weightValues.push_back(weightPtrs[1]);
154 strides.push_back(1);
155 strides.push_back(1);
156
157 try{
158 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
159 v_values, 1, true, weightPtrs[0], 1);
160 }
161 catch (std::exception &e){
162 fail = 1;
163 }
164
165 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
166
167 testBasisVector(ia, valueStrides, 1, 1);
168
169 delete ia;
170
171
172 }
173
174 {
175 // A Zoltan2::BasicVectorAdapter object with a multivector with weights
176
177 std::vector<const zscalar_t *> weightValues, values;
178 std::vector<int> wstrides, vstrides;
179
180 for (int dim=0; dim < wdim; dim++){
181 weightValues.push_back(weightPtrs[dim]);
182 wstrides.push_back(1);
183 }
184
185 for (int dim=0; dim < mvdim; dim++){
186 values.push_back(valuePtrs[dim]);
187 vstrides.push_back(mvdim);
188 }
189
190 try{
192 numLocalIds, myIds, values, vstrides, weightValues, wstrides);
193 }
194 catch (std::exception &e){
195 fail = 1;
196 }
197
198 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
199
200 testBasisVector(ia, valueStrides, wdim, mvdim);
201
202 delete ia;
203
204 }
205
206 if (rank == 0)
207 std::cout << "PASS" << std::endl;
208
209 delete [] myIds;
210 delete [] weights;
211 delete [] weightStrides;
212 delete [] weightPtrs;
213 delete [] v_values;
214 delete [] mv_values;
215 delete [] valueStrides;
216 delete [] valuePtrs;
217
218 return fail;
219}
220
void testBasisVector(Zoltan2::BasicVectorAdapter< userTypes_t > *ia, int *valueStrides, int wdim, int mvdim)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
Defines the BasicVectorAdapter class.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
int main()
Kokkos::View< scalar_t **, Kokkos::LayoutLeft, device_t > CoordsDeviceView
typename CoordsDeviceView::host_mirror_type CoordsHostView
Kokkos::View< const gno_t *, device_t > ConstIdsDeviceView
typename ConstIdsDeviceView::host_mirror_type ConstIdsHostView
typename WeightsDeviceView::host_mirror_type WeightsHostView
Kokkos::View< scalar_t **, device_t > WeightsDeviceView
A simple class that can be the User template argument for an InputAdapter.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
void getWeightsHostView(typename Base::WeightsHostView1D &hostWgts, int idx=0) const override
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getWeightsDeviceView(typename Base::WeightsDeviceView1D &deviceWgts, int idx=0) const override
void getIDsHostView(typename Base::ConstIdsHostView &ids) const override
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
void getIDsDeviceView(typename Base::ConstIdsDeviceView &ids) const override
void getIDsView(const gno_t *&ids) const
void getCoordinatesDeviceView(typename AdapterWithCoords< User >::CoordsDeviceView &elements) const override
void getCoordinatesView(const scalar_t *&elements, int &stride, int idx=0) const override
void getCoordinatesHostView(typename AdapterWithCoords< User >::CoordsHostView &elements) const override
static const std::string fail
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights