Zoltan2
Loading...
Searching...
No Matches
XpetraVectorInput.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
17#include <string>
18
22
23#include <Teuchos_DefaultComm.hpp>
24#include <Teuchos_RCP.hpp>
25#include <Teuchos_Comm.hpp>
26#include <Teuchos_CommHelpers.hpp>
27
28using Teuchos::RCP;
29using Teuchos::rcp;
30using Teuchos::rcp_const_cast;
31using Teuchos::Comm;
32
33typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> tvector_t;
34typedef Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> xvector_t;
35
36void printVector(RCP<const Comm<int> > &comm, zlno_t vlen,
37 const zgno_t *vtxIds, const zscalar_t *vals)
38{
39 int rank = comm->getRank();
40 int nprocs = comm->getSize();
41 comm->barrier();
42 for (int p=0; p < nprocs; p++){
43 if (p == rank){
44 std::cout << rank << ":" << std::endl;
45 for (zlno_t i=0; i < vlen; i++){
46 std::cout << " " << vtxIds[i] << ": " << vals[i] << std::endl;
47 }
48 std::cout.flush();
49 }
50 comm->barrier();
51 }
52 comm->barrier();
53}
54
55template <typename User>
58 zscalar_t **weights, int *strides)
59{
60 RCP<const Comm<int> > comm = vector.getMap()->getComm();
61 int fail = 0, gfail=0;
62
63 if (!fail && ia.getNumEntriesPerID() !=1)
64 fail = 42;
65
66 if (!fail && ia.getNumWeightsPerID() !=wdim)
67 fail = 41;
68
69 if (!fail && ia.getLocalNumIDs() != vector.getLocalLength())
70 fail = 4;
71
72 gfail = globalFail(*comm, fail);
73
74 if (!gfail){
75 const zgno_t *vtxIds=NULL;
76 const zscalar_t *vals=NULL;
77 int stride;
78
79 size_t nvals = ia.getLocalNumIDs();
80 if (nvals != vector.getLocalLength())
81 fail = 8;
82
83 ia.getIDsView(vtxIds);
84 ia.getEntriesView(vals, stride);
85
86 if (!fail && stride != 1)
87 fail = 10;
88
89 gfail = globalFail(*comm, fail);
90
91 if (gfail == 0){
92 printVector(comm, nvals, vtxIds, vals);
93 }
94 }
95
96 if (!gfail && wdim){
97 const zscalar_t *wgt =NULL;
98 int stride;
99
100 for (int w=0; !fail && w < wdim; w++){
101 ia.getWeightsView(wgt, stride, w);
102
103 if (!fail && stride != strides[w])
104 fail = 101;
105
106 for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
107 if (wgt[v*stride] != weights[w][v*stride])
108 fail=102;
109 }
110 }
111
112 gfail = globalFail(*comm, fail);
113 }
114
115 return gfail;
116}
117
118
119int main(int narg, char *arg[])
120{
121 Tpetra::ScopeGuard tscope(&narg, &arg);
122 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
123
124 int rank = comm->getRank();
125 int fail = 0, gfail=0;
126 bool aok = true;
127
128 // Create object that can give us test Tpetra, Xpetra
129 // vectors for testing.
130
131 RCP<UserInputForTests> uinput;
132 Teuchos::ParameterList params;
133 params.set("input file", "simple");
134 params.set("file type", "Chaco");
135 try{
136 uinput = rcp(new UserInputForTests(params, comm));
137 }
138 catch(std::exception &e){
139 aok = false;
140 std::cout << e.what() << std::endl;
141 }
142 TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
143
144 RCP<tvector_t> tV; // original vector (for checking)
145 RCP<tvector_t> newV; // migrated vector
146
147 tV = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
148 tV->randomize();
149 size_t vlen = tV->getLocalLength();
150
151 // To test migration in the input adapter we need a Solution object.
152
153 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
154
155 int nWeights = 1;
156
158 typedef adapter_t::part_t part_t;
159
160 part_t *p = new part_t [vlen];
161 memset(p, 0, sizeof(part_t) * vlen);
162 ArrayRCP<part_t> solnParts(p, 0, vlen, true);
163
164 std::vector<const zscalar_t *> emptyWeights;
165 std::vector<int> emptyStrides;
166
167 Zoltan2::PartitioningSolution<adapter_t> solution(env, comm, nWeights);
168 solution.setParts(solnParts);
169
171 // User object is Tpetra::Vector, no weights
172 if (!gfail){
173 if (rank==0)
174 std::cout << "Constructed with Tpetra::Vector" << std::endl;
175
176 RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
177 RCP<adapter_t> tVInput;
178
179 try {
180 tVInput =
182 emptyWeights, emptyStrides));
183 }
184 catch (std::exception &e){
185 aok = false;
186 std::cout << e.what() << std::endl;
187 }
188 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter ", 1);
189
190 fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
191
192 gfail = globalFail(*comm, fail);
193
194 if (!gfail){
195 tvector_t *vMigrate = NULL;
196 try{
197 tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
198 newV = rcp(vMigrate);
199 }
200 catch (std::exception &e){
201 fail = 11;
202 }
203
204 gfail = globalFail(*comm, fail);
205
206 if (!gfail){
207 RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
208 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
209 try{
210 newInput = rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(cnewV,
211 emptyWeights, emptyStrides));
212 }
213 catch (std::exception &e){
214 aok = false;
215 std::cout << e.what() << std::endl;
216 }
217 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 2 ", 1);
218
219 if (rank==0){
220 std::cout << "Constructed with ";
221 std::cout << "Tpetra::Vector migrated to proc 0" << std::endl;
222 }
223 fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
224 if (fail) fail += 100;
225 gfail = globalFail(*comm, fail);
226 }
227 }
228 if (gfail){
229 printFailureCode(*comm, fail);
230 }
231 }
232
234 // User object is Xpetra::Vector
235 if (!gfail){
236 if (rank==0)
237 std::cout << "Constructed with Xpetra::Vector" << std::endl;
238
239 RCP<tvector_t> xtV =
240 rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
241 xtV->randomize();
243 RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
244 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
245
246 try {
247 xVInput =
249 emptyWeights, emptyStrides));
250 }
251 catch (std::exception &e){
252 aok = false;
253 std::cout << e.what() << std::endl;
254 }
255 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 3 ", 1);
256
257 fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
258
259 gfail = globalFail(*comm, fail);
260
261 if (!gfail){
262 xvector_t *vMigrate =NULL;
263 try{
264 xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
265 }
266 catch (std::exception &e){
267 fail = 11;
268 }
269
270 gfail = globalFail(*comm, fail);
271
272 if (!gfail){
273 RCP<const xvector_t> cnewV(vMigrate);
274 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
275 try{
276 newInput =
278 emptyWeights, emptyStrides));
279 }
280 catch (std::exception &e){
281 aok = false;
282 std::cout << e.what() << std::endl;
283 }
284 TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 4 ", 1);
285
286 if (rank==0){
287 std::cout << "Constructed with ";
288 std::cout << "Xpetra::Vector migrated to proc 0" << std::endl;
289 }
290 fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
291 if (fail) fail += 100;
292 gfail = globalFail(*comm, fail);
293 }
294 }
295 if (gfail){
296 printFailureCode(*comm, fail);
297 }
298 }
299
301 // DONE
302
303 if (rank==0)
304 std::cout << "PASS" << std::endl;
305}
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
Xpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Tpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int wdim, zscalar_t **weights, int *strides)
void printVector(RCP< const Comm< int > > &comm, zlno_t vlen, const zgno_t *vtxIds, const zscalar_t *vals)
Traits for application input objects.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
Defines the XpetraMultiVectorAdapter.
int main()
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
A PartitioningSolution is a solution to a partitioning problem.
void setParts(ArrayRCP< part_t > &partList)
The algorithm uses setParts to set the solution.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
int getNumEntriesPerID() const
Return the number of vectors.
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
size_t getLocalNumIDs() const
Returns the number of objects on this process.
static const std::string fail
static ArrayRCP< ArrayRCP< zscalar_t > > weights
SparseMatrixAdapter_t::part_t part_t
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.