50 #include <Teuchos_RCP.hpp> 51 #include <Tpetra_CrsMatrix.hpp> 52 #include <Tpetra_DefaultPlatform.hpp> 53 #include <MatrixMarket_Tpetra.hpp> 63 typedef Tpetra::Map<zlno_t, zgno_t>
zmap_t;
64 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t>
zmatrix_t;
65 typedef Tpetra::CrsGraph<zlno_t, zgno_t>
zgraph_t;
71 template <
typename CC_T>
83 if (cc.getNumComponents() != nccAnswer) {
84 std::cout << name <<
"Invalid number of components " 85 << cc.getNumComponents() <<
" should be " << nccAnswer
90 if (cc.getMaxComponentSize() != maxAnswer) {
91 std::cout << name <<
"Maximum component size " 92 << cc.getMaxComponentSize() <<
" should be " << maxAnswer
97 if (cc.getMinComponentSize() != minAnswer) {
98 std::cout << name <<
"Minimum component size " 99 << cc.getMinComponentSize() <<
" should be " << minAnswer
104 if (cc.getAvgComponentSize() != avgAnswer) {
105 std::cout << name <<
"Average component size " 106 << cc.getAvgComponentSize() <<
" should be " << avgAnswer
120 std::ostringstream namestream;
121 namestream << comm->getRank() <<
" every_third ";
122 std::string name = namestream.str();
127 const size_t gNvtx = 27;
129 Teuchos::RCP<const zmap_t> map = rcp(
new zmap_t(gNvtx, 0, comm));
130 size_t nVtx = map->getNodeNumElements();
133 size_t maxRowLen = 1;
134 Teuchos::RCP<zmatrix_t> matrix = rcp(
new zmatrix_t(map, maxRowLen));
136 Teuchos::Array<zgno_t> col(3);
137 Teuchos::Array<zscalar_t> val(3); val[0] = 1.; val[1] = 1.; val[2] = 1.;
139 for (
size_t i = 0; i < nVtx; i++) {
140 zgno_t id = map->getGlobalElement(i);
141 col[0] = (
id+3)%gNvtx;
142 col[1] = (
id+6)%gNvtx;
143 col[2] = (
id+9)%gNvtx;
144 matrix->insertGlobalValues(
id, col(), val());
147 matrix->fillComplete(map, map);
157 size_t nccAnswer = std::min<size_t>(3, nVtx);
158 size_t maxAnswer = nVtx/3 + ((nVtx%3)!=0);
159 size_t minAnswer = (nVtx ? std::max<size_t>(nVtx/3, 1) : 0);
160 double avgAnswer = (nVtx ? double(nVtx) / double(std::min<size_t>(nVtx,3))
163 ierr = checkResult<cc_t>(name, cc,
164 nccAnswer, maxAnswer, minAnswer, avgAnswer);
175 std::ostringstream namestream;
176 namestream << comm->getRank() <<
" dist_component ";
177 std::string name = namestream.str();
182 const size_t gNvtx = 25;
184 Teuchos::RCP<const zmap_t> map = rcp(
new zmap_t(gNvtx, 0, comm));
185 size_t nVtx = map->getNodeNumElements();
188 size_t maxRowLen = 1;
189 Teuchos::RCP<zmatrix_t> matrix = rcp(
new zmatrix_t(map, maxRowLen));
191 Teuchos::Array<zgno_t> col(3);
192 Teuchos::Array<zscalar_t> val(3); val[0] = 1.; val[1] = 1.; val[2] = 1.;
194 for (
size_t i = 0; i < nVtx; i++) {
195 zgno_t id = map->getGlobalElement(i);
196 col[0] = (
id+4)%gNvtx;
197 col[1] = (
id+1)%gNvtx;
198 col[2] = (
id+3)%gNvtx;
199 matrix->insertGlobalValues(
id, col(), val());
202 matrix->fillComplete(map, map);
213 size_t nccAnswer = size_t(nVtx > 0);
214 size_t maxAnswer = nVtx;
215 size_t minAnswer = nVtx;
216 double avgAnswer = nVtx;
218 ierr = checkResult<cc_t>(name, cc,
219 nccAnswer, maxAnswer, minAnswer, avgAnswer);
230 std::ostringstream namestream;
231 namestream << comm->getRank() <<
" one_proc ";
232 std::string name = namestream.str();
236 int me = comm->getRank();
237 int np = comm->getSize();
238 bool allOnThisProc = (me == np-1);
241 const size_t gNvtx = 25;
244 if (allOnThisProc) nVtx = gNvtx;
247 Teuchos::RCP<const zmap_t> map = rcp(
new zmap_t(gNvtx, nVtx, 0, comm));
250 size_t maxRowLen = 1;
251 Teuchos::RCP<zmatrix_t> matrix = rcp(
new zmatrix_t(map, maxRowLen));
253 Teuchos::Array<zgno_t> col(2);
254 Teuchos::Array<zscalar_t> val(2); val[0] = 1.; val[1] = 1.;
255 for (
size_t i = 0; i < nVtx; i++) {
256 zgno_t id = map->getGlobalElement(i);
258 col[1] = (
id+1)%gNvtx;
259 matrix->insertGlobalValues(
id, col(), val());
262 matrix->fillComplete(map, map);
274 size_t nccAnswer = 1;
275 size_t maxAnswer = gNvtx;
276 size_t minAnswer = gNvtx;
277 double avgAnswer = double(gNvtx);
279 ierr = checkResult<cc_t>(name, cc,
280 nccAnswer, maxAnswer, minAnswer, avgAnswer);
284 size_t nccAnswer = 0;
285 size_t maxAnswer = 0;
286 size_t minAnswer = 0;
287 double avgAnswer = 0.;
289 ierr = checkResult<cc_t>(name, cc,
290 nccAnswer, maxAnswer, minAnswer, avgAnswer);
303 std::ostringstream namestream;
304 namestream << comm->getRank() <<
" no_graph ";
305 std::string name = namestream.str();
310 const size_t gNvtx = 25;
312 Teuchos::RCP<const zmap_t> map = rcp(
new zmap_t(gNvtx, 0, comm));
313 size_t nVtx = map->getNodeNumElements();
316 size_t maxRowLen = 1;
317 Teuchos::RCP<zmatrix_t> matrix = rcp(
new zmatrix_t(map, maxRowLen));
318 matrix->fillComplete(map, map);
323 #ifdef HAVE_ZOLTAN2_MPI 327 if (comm->getRank() == 0) std::cout <<
" using MPI_Comm " << std::endl;
328 MPI_Comm useThisComm = Teuchos::getRawMpiComm(*comm);
330 const Teuchos::Comm<int> &useThisComm = *comm;
335 cc_t cc(ia, useThisComm);
339 size_t nccAnswer = nVtx;
340 size_t maxAnswer = size_t(nVtx > 0);
341 size_t minAnswer = size_t(nVtx > 0);
342 double avgAnswer = double(nVtx > 0);
344 ierr = checkResult<cc_t>(name, cc,
345 nccAnswer, maxAnswer, minAnswer, avgAnswer);
354 Teuchos::GlobalMPISession mpiSession(&narg, &arg, NULL);
355 Teuchos::RCP<const Teuchos::Comm<int> > comm =
356 Tpetra::DefaultPlatform::getDefaultPlatform().getComm();
357 int me = comm->getRank();
360 if (me == 0) std::cout <<
"test_one_proc..." << std::endl;
362 if (me == 0) std::cout <<
"test_no_graph..." << std::endl;
364 if (me == 0) std::cout <<
"test_dist_component..." << std::endl;
366 if (me == 0) std::cout <<
"test_every_third..." << std::endl;
370 Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1,
371 &testReturn, >estReturn);
373 if (gtestReturn) std::cout <<
"FAIL" << std::endl;
374 else std::cout <<
"PASS" << std::endl;
Tpetra::CrsGraph< zlno_t, zgno_t > zgraph_t
Provides access for Zoltan2 to Xpetra::CrsMatrix data.
int test_dist_component(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Zoltan2::XpetraCrsMatrixAdapter< zmatrix_t > matrixAdapter_t
int checkResult(std::string &name, CC_T cc, size_t nccAnswer, size_t maxAnswer, size_t minAnswer, double avgAnswer)
Provides access for Zoltan2 to Xpetra::CrsGraph data.
Zoltan2::XpetraCrsGraphAdapter< zgraph_t > graphAdapter_t
int test_every_third(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
common code used by tests
int test_one_proc(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Tpetra::CrsMatrix< zscalar_t, zlno_t, zgno_t > zmatrix_t
Defines XpetraCrsGraphAdapter class.
Defines the XpetraCrsMatrixAdapter class.
int test_no_graph(Teuchos::RCP< const Teuchos::Comm< int > > &comm)
int main(int narg, char **arg)
Identify and compute the number of connected components in a processor's input Note that this routine...
Tpetra::Map< zlno_t, zgno_t > zmap_t