Zoltan2
Toggle main menu visibility
Loading...
Searching...
No Matches
componentMetrics.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
#include <
Zoltan2_componentMetrics.hpp
>
46
#include <
Zoltan2_XpetraCrsMatrixAdapter.hpp
>
47
#include <
Zoltan2_XpetraCrsGraphAdapter.hpp
>
48
#include <
Zoltan2_TestHelpers.hpp
>
49
#include <iostream>
50
#include <Teuchos_RCP.hpp>
51
#include <Tpetra_CrsMatrix.hpp>
52
#include <MatrixMarket_Tpetra.hpp>
53
54
56
// Test program for componentMetrics code.
57
// Usage:
58
// a.out
59
// Karen Devine, 2017
61
62
typedef
Tpetra::Map<zlno_t, zgno_t>
zmap_t
;
63
typedef
Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t>
zmatrix_t
;
64
typedef
Tpetra::CrsGraph<zlno_t, zgno_t>
zgraph_t
;
65
66
typedef
Zoltan2::XpetraCrsMatrixAdapter<zmatrix_t>
matrixAdapter_t
;
67
typedef
Zoltan2::XpetraCrsGraphAdapter<zgraph_t>
graphAdapter_t
;
68
70
template
<
typename
CC_T>
71
int
checkResult
(
72
std::string &name,
// test name (including rank)
73
CC_T cc,
// the perProcessorComponentMetrics object for the test
74
size_t
nccAnswer,
// Expected answers for the test
75
size_t
maxAnswer,
76
size_t
minAnswer,
77
double
avgAnswer
78
)
79
{
80
int
ierr = 0;
81
82
if
(cc.getNumComponents() != nccAnswer) {
83
std::cout << name <<
"Invalid number of components "
84
<< cc.getNumComponents() <<
" should be "
<< nccAnswer
85
<< std::endl;
86
ierr++;
87
}
88
89
if
(cc.getMaxComponentSize() != maxAnswer) {
90
std::cout << name <<
"Maximum component size "
91
<< cc.getMaxComponentSize() <<
" should be "
<< maxAnswer
92
<< std::endl;
93
ierr++;
94
}
95
96
if
(cc.getMinComponentSize() != minAnswer) {
97
std::cout << name <<
"Minimum component size "
98
<< cc.getMinComponentSize() <<
" should be "
<< minAnswer
99
<< std::endl;
100
ierr++;
101
}
102
103
if
(cc.getAvgComponentSize() != avgAnswer) {
104
std::cout << name <<
"Average component size "
105
<< cc.getAvgComponentSize() <<
" should be "
<< avgAnswer
106
<< std::endl;
107
ierr++;
108
}
109
110
return
ierr;
111
}
112
114
int
test_every_third
(Teuchos::RCP<
const
Teuchos::Comm<int> > &comm)
115
{
116
// Testing a graph with every third vertex in the same component
117
// Test using a GraphAdapter.
118
119
std::ostringstream namestream;
120
namestream << comm->getRank() <<
" every_third "
;
121
std::string name = namestream.str();
122
123
int
ierr = 0;
124
125
// Create a default map
126
const
size_t
gNvtx = 27;
127
128
Teuchos::RCP<const zmap_t> map = rcp(
new
zmap_t
(gNvtx, 0, comm));
129
size_t
nVtx = map->getLocalNumElements();
130
131
// Create a Tpetra::Matrix with every third local row in the same component
132
size_t
maxRowLen = 3;
133
Teuchos::RCP<zmatrix_t> matrix = rcp(
new
zmatrix_t
(map, maxRowLen));
134
135
Teuchos::Array<zgno_t> col(3);
136
Teuchos::Array<zscalar_t> val(3); val[0] = 1.; val[1] = 1.; val[2] = 1.;
137
138
for
(
size_t
i = 0; i < nVtx; i++) {
139
zgno_t
id
= map->getGlobalElement(i);
140
col[0] = (
id
+3)%gNvtx;
141
col[1] = (
id
+6)%gNvtx;
142
col[2] = (
id
+9)%gNvtx;
143
matrix->insertGlobalValues(
id
, col(), val());
144
}
145
146
matrix->fillComplete(map, map);
147
148
// Create a Zoltan2 XpetraGraphAdapter
149
graphAdapter_t
ia(matrix->getCrsGraph(), 0);
150
151
// Call connected components utility
152
typedef
Zoltan2::perProcessorComponentMetrics<graphAdapter_t>
cc_t;
153
cc_t cc(ia, *comm);
154
155
// Check result:
156
size_t
nccAnswer = std::min<size_t>(3, nVtx);
157
size_t
maxAnswer = nVtx/3 + ((nVtx%3)!=0);
158
size_t
minAnswer = (nVtx ? std::max<size_t>(nVtx/3, 1) : 0);
159
double
avgAnswer = (nVtx ? double(nVtx) / double(std::min<size_t>(nVtx,3))
160
: 0.);
161
162
ierr =
checkResult<cc_t>
(name, cc,
163
nccAnswer, maxAnswer, minAnswer, avgAnswer);
164
165
return
ierr;
166
}
167
169
int
test_dist_component
(Teuchos::RCP<
const
Teuchos::Comm<int> > &comm)
170
{
171
// Testing a graph with single component, distributed among procs
172
// Test using a GraphAdapter.
173
174
std::ostringstream namestream;
175
namestream << comm->getRank() <<
" dist_component "
;
176
std::string name = namestream.str();
177
178
int
ierr = 0;
179
180
// Create a default map
181
const
size_t
gNvtx = 25;
182
183
Teuchos::RCP<const zmap_t> map = rcp(
new
zmap_t
(gNvtx, 0, comm));
184
size_t
nVtx = map->getLocalNumElements();
185
186
// Create a Tpetra::Matrix with a single component
187
size_t
maxRowLen = 3;
188
Teuchos::RCP<zmatrix_t> matrix = rcp(
new
zmatrix_t
(map, maxRowLen));
189
190
Teuchos::Array<zgno_t> col(3);
191
Teuchos::Array<zscalar_t> val(3); val[0] = 1.; val[1] = 1.; val[2] = 1.;
192
193
for
(
size_t
i = 0; i < nVtx; i++) {
194
zgno_t
id
= map->getGlobalElement(i);
195
col[0] = (
id
+4)%gNvtx;
196
col[1] = (
id
+1)%gNvtx;
197
col[2] = (
id
+3)%gNvtx;
198
matrix->insertGlobalValues(
id
, col(), val());
199
}
200
201
matrix->fillComplete(map, map);
202
203
// Create a Zoltan2 XpetraGraphAdapter
204
graphAdapter_t
ia(matrix->getCrsGraph(), 0);
205
206
// Call connected components utility
207
typedef
Zoltan2::perProcessorComponentMetrics<graphAdapter_t>
cc_t;
208
cc_t cc(ia, *comm);
209
210
// Check result:
211
// one component with all local vertices
212
size_t
nccAnswer = size_t(nVtx > 0);
213
size_t
maxAnswer = nVtx;
214
size_t
minAnswer = nVtx;
215
double
avgAnswer = nVtx;
216
217
ierr =
checkResult<cc_t>
(name, cc,
218
nccAnswer, maxAnswer, minAnswer, avgAnswer);
219
220
return
ierr;
221
}
222
224
int
test_one_proc
(Teuchos::RCP<
const
Teuchos::Comm<int> > &comm)
225
{
226
// Testing a graph with all vertices and edges on a single proc
227
// Test using a MatrixAdapter.
228
229
std::ostringstream namestream;
230
namestream << comm->getRank() <<
" one_proc "
;
231
std::string name = namestream.str();
232
233
int
ierr = 0;
234
235
int
me = comm->getRank();
236
int
np = comm->getSize();
237
bool
allOnThisProc = (me == np-1);
238
239
// Create a map with all data on a single proc
240
const
size_t
gNvtx = 25;
241
size_t
nVtx;
242
243
if
(allOnThisProc) nVtx = gNvtx;
244
else
nVtx = 0;
245
246
Teuchos::RCP<const zmap_t> map = rcp(
new
zmap_t
(gNvtx, nVtx, 0, comm));
247
248
// Create a Tpetra::Matrix with one component
249
size_t
maxRowLen = 2;
250
Teuchos::RCP<zmatrix_t> matrix = rcp(
new
zmatrix_t
(map, maxRowLen));
251
if
(allOnThisProc) {
252
Teuchos::Array<zgno_t> col(2);
253
Teuchos::Array<zscalar_t> val(2); val[0] = 1.; val[1] = 1.;
254
for
(
size_t
i = 0; i < nVtx; i++) {
255
zgno_t
id
= map->getGlobalElement(i);
256
col[0] = id;
257
col[1] = (
id
+1)%gNvtx;
258
matrix->insertGlobalValues(
id
, col(), val());
259
}
260
}
261
matrix->fillComplete(map, map);
262
263
// Create a Zoltan2 XpetraMatrixAdapter
264
matrixAdapter_t
ia(matrix, 0);
265
266
// Call connected components utility
267
typedef
Zoltan2::perProcessorComponentMetrics<matrixAdapter_t>
cc_t;
268
cc_t cc(ia, *comm);
269
270
// Check result:
271
if
(allOnThisProc) {
272
// one component with all global vertices
273
size_t
nccAnswer = 1;
274
size_t
maxAnswer = gNvtx;
275
size_t
minAnswer = gNvtx;
276
double
avgAnswer = double(gNvtx);
277
278
ierr =
checkResult<cc_t>
(name, cc,
279
nccAnswer, maxAnswer, minAnswer, avgAnswer);
280
}
281
else
{
282
// one component with all global vertices
283
size_t
nccAnswer = 0;
284
size_t
maxAnswer = 0;
285
size_t
minAnswer = 0;
286
double
avgAnswer = 0.;
287
288
ierr =
checkResult<cc_t>
(name, cc,
289
nccAnswer, maxAnswer, minAnswer, avgAnswer);
290
}
291
292
return
ierr;
293
}
294
296
int
test_no_graph
(Teuchos::RCP<
const
Teuchos::Comm<int> > &comm)
297
{
298
// Testing a graph with no edges
299
// Test using a MatrixAdapter.
300
// Test using an MPI communicator rather than Teuchos::Comm, if appropriate.
301
302
std::ostringstream namestream;
303
namestream << comm->getRank() <<
" no_graph "
;
304
std::string name = namestream.str();
305
306
int
ierr = 0;
307
308
// Create a default map
309
const
size_t
gNvtx = 25;
310
311
Teuchos::RCP<const zmap_t> map = rcp(
new
zmap_t
(gNvtx, 0, comm));
312
size_t
nVtx = map->getLocalNumElements();
313
314
// Create a Tpetra::Matrix with no edges
315
size_t
maxRowLen = 1;
316
Teuchos::RCP<zmatrix_t> matrix = rcp(
new
zmatrix_t
(map, maxRowLen));
317
matrix->fillComplete(map, map);
318
319
// Create a Zoltan2 XpetraMatrixAdapter
320
matrixAdapter_t
ia(matrix, 0);
321
322
#ifdef HAVE_ZOLTAN2_MPI
323
// For testing only, extract MPI communicator from Teuchos::Comm
324
// There's no need to do this in real life; use Teuchos::Comm if you have it.
325
// We just want to exercise the MPI_Comm interface here.
326
if
(comm->getRank() == 0) std::cout <<
" using MPI_Comm "
<< std::endl;
327
MPI_Comm useThisComm = Teuchos::getRawMpiComm(*comm);
328
#else
329
const
Teuchos::Comm<int> &useThisComm = *comm;
330
#endif
331
332
// Call connected components utility
333
typedef
Zoltan2::perProcessorComponentMetrics<matrixAdapter_t>
cc_t;
334
cc_t cc(ia, useThisComm);
335
336
// Check result:
337
// With no edges, every vertex should be a component
338
size_t
nccAnswer = nVtx;
339
size_t
maxAnswer = size_t(nVtx > 0);
340
size_t
minAnswer = size_t(nVtx > 0);
341
double
avgAnswer = double(nVtx > 0);
342
343
ierr =
checkResult<cc_t>
(name, cc,
344
nccAnswer, maxAnswer, minAnswer, avgAnswer);
345
346
return
ierr;
347
}
348
350
int
main
(
int
narg,
char
** arg)
351
{
352
// Establish session.
353
Tpetra::ScopeGuard tscope(&narg, &arg);
354
Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
355
int
me = comm->getRank();
356
357
int
testReturn = 0;
358
if
(me == 0) std::cout <<
"test_one_proc..."
<< std::endl;
359
testReturn +=
test_one_proc
(comm);
// all data on a single proc
360
if
(me == 0) std::cout <<
"test_no_graph..."
<< std::endl;
361
testReturn +=
test_no_graph
(comm);
// no edges in graph
362
if
(me == 0) std::cout <<
"test_dist_component..."
<< std::endl;
363
testReturn +=
test_dist_component
(comm);
// one component per rank
364
if
(me == 0) std::cout <<
"test_every_third..."
<< std::endl;
365
testReturn +=
test_every_third
(comm);
// every 3rd vtx in same component
366
367
int
gtestReturn = 0;
368
Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1,
369
&testReturn, >estReturn);
370
if
(me == 0) {
371
if
(gtestReturn) std::cout <<
"FAIL"
<< std::endl;
372
else
std::cout <<
"PASS"
<< std::endl;
373
}
374
375
return
gtestReturn;
376
}
Zoltan2_TestHelpers.hpp
common code used by tests
zgno_t
Tpetra::Map ::global_ordinal_type zgno_t
Definition
Zoltan2_TestHelpers.hpp:107
Zoltan2_XpetraCrsGraphAdapter.hpp
Defines XpetraCrsGraphAdapter class.
Zoltan2_XpetraCrsMatrixAdapter.hpp
Defines the XpetraCrsMatrixAdapter class.
Zoltan2_componentMetrics.hpp
Identify and compute the number of connected components in a processor's input Note that this routine...
main
int main()
Definition
absdefinitiontest.cpp:6
Zoltan2::XpetraCrsGraphAdapter
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Definition
Zoltan2_XpetraCrsGraphAdapter.hpp:84
Zoltan2::XpetraCrsMatrixAdapter
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
Definition
Zoltan2_XpetraCrsMatrixAdapter.hpp:86
Zoltan2::perProcessorComponentMetrics
Definition
Zoltan2_componentMetrics.hpp:65
zmap_t
Tpetra::Map< zlno_t, zgno_t > zmap_t
Definition
componentMetrics.cpp:62
zgraph_t
Tpetra::CrsGraph< zlno_t, zgno_t > zgraph_t
Definition
componentMetrics.cpp:64
test_dist_component
int test_dist_component(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Definition
componentMetrics.cpp:169
test_one_proc
int test_one_proc(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Definition
componentMetrics.cpp:224
test_every_third
int test_every_third(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Definition
componentMetrics.cpp:114
test_no_graph
int test_no_graph(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Definition
componentMetrics.cpp:296
graphAdapter_t
Zoltan2::XpetraCrsGraphAdapter< zgraph_t > graphAdapter_t
Definition
componentMetrics.cpp:67
checkResult
int checkResult(std::string &name, CC_T cc, size_t nccAnswer, size_t maxAnswer, size_t minAnswer, double avgAnswer)
Definition
componentMetrics.cpp:71
zmatrix_t
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t > zmatrix_t
Definition
componentMetrics.cpp:63
matrixAdapter_t
Zoltan2::XpetraCrsMatrixAdapter< tMatrix_t, tMVector_t > matrixAdapter_t
Definition
zoltanCompare.cpp:120
test
core
unit
util
componentMetrics.cpp
Generated by
1.17.0