|
Zoltan2
|
00001 // @HEADER 00002 // 00003 // *********************************************************************** 00004 // 00005 // Zoltan2: A package of combinatorial algorithms for scientific computing 00006 // Copyright 2012 Sandia Corporation 00007 // 00008 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation, 00009 // the U.S. Government retains certain rights in this software. 00010 // 00011 // Redistribution and use in source and binary forms, with or without 00012 // modification, are permitted provided that the following conditions are 00013 // met: 00014 // 00015 // 1. Redistributions of source code must retain the above copyright 00016 // notice, this list of conditions and the following disclaimer. 00017 // 00018 // 2. Redistributions in binary form must reproduce the above copyright 00019 // notice, this list of conditions and the following disclaimer in the 00020 // documentation and/or other materials provided with the distribution. 00021 // 00022 // 3. Neither the name of the Corporation nor the names of the 00023 // contributors may be used to endorse or promote products derived from 00024 // this software without specific prior written permission. 00025 // 00026 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY 00027 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE 00028 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR 00029 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE 00030 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, 00031 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, 00032 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR 00033 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF 00034 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING 00035 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS 00036 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE. 00037 // 00038 // Questions? Contact Karen Devine (kddevin@sandia.gov) 00039 // Erik Boman (egboman@sandia.gov) 00040 // Siva Rajamanickam (srajama@sandia.gov) 00041 // 00042 // *********************************************************************** 00043 // 00044 // @HEADER 00045 // 00046 // Testing of GraphModel built from Xpetra matrix input adapters. 00047 // 00048 00058 #include <Zoltan2_GraphModel.hpp> 00059 #include <Zoltan2_XpetraCrsMatrixAdapter.hpp> 00060 #include <Zoltan2_XpetraCrsGraphAdapter.hpp> 00061 #include <Zoltan2_BasicVectorAdapter.hpp> 00062 #include <Zoltan2_TestHelpers.hpp> 00063 #include <Zoltan2_InputTraits.hpp> 00064 00065 #include <string> 00066 #include <vector> 00067 #include <iostream> 00068 #include <bitset> 00069 00070 #include <Teuchos_Comm.hpp> 00071 #include <Teuchos_DefaultComm.hpp> 00072 #include <Teuchos_ArrayView.hpp> 00073 00074 using namespace std; 00075 using Teuchos::RCP; 00076 using Teuchos::rcp; 00077 using Teuchos::Comm; 00078 using Teuchos::DefaultComm; 00079 using Teuchos::ArrayView; 00080 00081 typedef Zoltan2::BasicUserTypes<zscalar_t, zzgid_t, zlno_t, zgno_t> simpleUser_t; 00082 00083 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> tcrsMatrix_t; 00084 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> tcrsGraph_t; 00085 typedef Tpetra::Map<zlno_t, zgno_t, znode_t> tmap_t; 00086 00087 typedef Zoltan2::StridedData<zlno_t, zscalar_t> input_t; 00088 00089 typedef Zoltan2::BasicVectorAdapter<simpleUser_t> simpleVAdapter_t; 00090 00091 typedef Zoltan2::MatrixAdapter<tcrsMatrix_t,simpleUser_t> baseMAdapter_t; 00092 typedef Zoltan2::GraphAdapter<tcrsGraph_t,simpleUser_t> baseGAdapter_t; 00093 00094 typedef Zoltan2::XpetraCrsMatrixAdapter<tcrsMatrix_t,simpleUser_t> xMAdapter_t; 00095 typedef Zoltan2::XpetraCrsGraphAdapter<tcrsGraph_t,simpleUser_t> xGAdapter_t; 00096 00097 using std::string; 00098 using std::vector; 00099 00101 void printGraph(zlno_t nrows, const zgno_t *v, 00102 const zlno_t *elid, const zgno_t *egid, 00103 const int *owner, const zlno_t *idx, 00104 const RCP<const Comm<int> > &comm) 00105 { 00106 int rank = comm->getRank(); 00107 int nprocs = comm->getSize(); 00108 comm->barrier(); 00109 if (rank == 0){ 00110 if (owner) 00111 std::cout << "Global graph:" << std::endl; 00112 else 00113 std::cout << "Local graph:" << std::endl; 00114 } 00115 00116 for (int p=0; p < nprocs; p++){ 00117 if (p == rank){ 00118 std::cout << "Rank " << p << std::endl; 00119 if (owner){ 00120 for (zlno_t i=0; i < nrows; i++){ 00121 std::cout << " Vtx " << *v++ << ": "; 00122 if (elid) 00123 for (zlno_t j=idx[i]; j < idx[i+1]; j++) 00124 std::cout << *elid++ << " (" << *owner++ << ") "; 00125 else 00126 for (zlno_t j=idx[i]; j < idx[i+1]; j++) 00127 std::cout << *egid++ << " (" << *owner++ << ") "; 00128 std::cout << std::endl; 00129 } 00130 std::cout.flush(); 00131 } 00132 else{ 00133 for (zlno_t i=0; i < nrows; i++){ 00134 std::cout << " Vtx " << i << ": "; 00135 if (elid) 00136 for (zlno_t j=idx[i]; j < idx[i+1]; j++) 00137 std::cout << *elid++ << " "; 00138 else 00139 for (zlno_t j=idx[i]; j < idx[i+1]; j++) 00140 std::cout << *egid++ << " "; 00141 std::cout << std::endl; 00142 } 00143 std::cout.flush(); 00144 } 00145 } 00146 comm->barrier(); 00147 } 00148 comm->barrier(); 00149 } 00150 00152 template <typename BaseAdapter, typename Adapter, typename MatrixOrGraph> 00153 void testAdapter( 00154 RCP<const MatrixOrGraph> &M, 00155 RCP<const Tpetra::CrsGraph<zlno_t, zgno_t> > &Mgraph, 00156 const RCP<const Comm<int> > &comm, 00157 bool idsAreConsecutive, 00158 int nVtxWeights, int nEdgeWeights, int nnzWgtIdx, int coordDim, 00159 bool consecutiveIdsRequested, bool removeSelfEdges) 00160 { 00161 int fail=0; 00162 int rank = comm->getRank(); 00163 int nprocs = comm->getSize(); 00164 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment); 00165 00166 zlno_t nLocalRows = M->getNodeNumRows(); 00167 zlno_t nLocalNZ = M->getNodeNumEntries(); 00168 zgno_t nGlobalRows = M->getGlobalNumRows(); 00169 zgno_t nGlobalNZ = M->getGlobalNumEntries(); 00170 00171 std::bitset<Zoltan2::NUM_MODEL_FLAGS> modelFlags; 00172 if (consecutiveIdsRequested) 00173 modelFlags.set(Zoltan2::IDS_MUST_BE_GLOBALLY_CONSECUTIVE); 00174 if (removeSelfEdges) 00175 modelFlags.set(Zoltan2::SELF_EDGES_MUST_BE_REMOVED); 00176 00177 // Set up some fake input 00178 zscalar_t **coords=NULL; 00179 zscalar_t **rowWeights=NULL; 00180 00181 if (coordDim > 0){ 00182 coords = new zscalar_t * [coordDim]; 00183 for (int i=0; i < coordDim; i++){ 00184 coords[i] = new zscalar_t [nLocalRows]; 00185 for (zlno_t j=0; j < nLocalRows; j++){ 00186 coords[i][j] = 100000*i + j; 00187 } 00188 } 00189 } 00190 00191 if (nVtxWeights > 0){ 00192 rowWeights = new zscalar_t * [nVtxWeights]; 00193 for (int i=0; i < nVtxWeights; i++){ 00194 if (nnzWgtIdx == i) 00195 rowWeights[i] = NULL; 00196 else{ 00197 rowWeights[i] = new zscalar_t [nLocalRows]; 00198 for (zlno_t j=0; j < nLocalRows; j++){ 00199 rowWeights[i][j] = 200000*i + j; 00200 } 00201 } 00202 } 00203 } 00204 00205 if (nEdgeWeights > 0){ 00206 printf("TODO: STILL NEED TO TEST EDGE WEIGHTS!\n"); 00207 } 00208 00209 // Create a matrix or graph input adapter. 00210 00211 Adapter tmi(M, nVtxWeights); 00212 00213 for (int i=0; i < nVtxWeights; i++){ 00214 if (nnzWgtIdx == i) 00215 tmi.setWeightIsDegree(i); 00216 else 00217 tmi.setWeights(rowWeights[i], 1, i); 00218 } 00219 00220 zzgid_t *gids = NULL; 00221 00222 simpleVAdapter_t *via = NULL; 00223 00224 if (coordDim > 0) { 00225 gids = new zzgid_t[nLocalRows]; 00226 for (zlno_t i = 0; i < nLocalRows; i++) 00227 gids[i] = M->getRowMap()->getGlobalElement(i); 00228 via = new simpleVAdapter_t(nLocalRows, gids, coords[0], 00229 (coordDim > 1 ? coords[1] : NULL), 00230 (coordDim > 2 ? coords[2] : NULL), 00231 1,1,1); 00232 tmi.setCoordinateInput(via); 00233 } 00234 00235 int numLocalDiags = M->getNodeNumDiags(); 00236 int numGlobalDiags = M->getGlobalNumDiags(); 00237 00238 const RCP<const tmap_t> rowMap = M->getRowMap(); 00239 const RCP<const tmap_t> colMap = M->getColMap(); 00240 00241 // How many neighbor vertices are not on my process. 00242 00243 int *numNbors = new int [nLocalRows]; 00244 int *numLocalNbors = new int [nLocalRows]; 00245 bool *haveDiag = new bool [nLocalRows]; 00246 zgno_t totalLocalNbors = 0; 00247 00248 for (zlno_t i=0; i < nLocalRows; i++){ 00249 numLocalNbors[i] = 0; 00250 haveDiag[i] = false; 00251 ArrayView<const zlno_t> idx; 00252 Mgraph->getLocalRowView(i, idx); 00253 numNbors[i] = idx.size(); 00254 00255 for (zlno_t j=0; j < idx.size(); j++){ 00256 if (idx[j] == i){ 00257 haveDiag[i] = true; 00258 numLocalNbors[i]++; 00259 totalLocalNbors++; 00260 } 00261 else{ 00262 zgno_t gidVal = colMap->getGlobalElement(idx[j]); 00263 if (rowMap->getLocalElement(gidVal) != 00264 Teuchos::OrdinalTraits<zlno_t>::invalid()) { 00265 // nbor is local to this process 00266 numLocalNbors[i]++; 00267 totalLocalNbors++; 00268 } 00269 } 00270 } 00271 } 00272 00273 // Create a GraphModel based on this input data. 00274 00275 if (rank == 0) std::cout << " Creating GraphModel" << std::endl; 00276 Zoltan2::GraphModel<BaseAdapter> *model = NULL; 00277 const BaseAdapter *baseTmi = dynamic_cast<BaseAdapter *>(&tmi); 00278 00279 try{ 00280 model = new Zoltan2::GraphModel<BaseAdapter>(baseTmi, env, 00281 comm, modelFlags); 00282 } 00283 catch (std::exception &e){ 00284 std::cerr << rank << ") " << e.what() << std::endl; 00285 fail = 1; 00286 } 00287 TEST_FAIL_AND_EXIT(*comm, !fail, "Creating graph model", 1) 00288 00289 // Test the GraphModel interface 00290 00291 if (rank == 0) std::cout << " Checking counts" << std::endl; 00292 if (model->getLocalNumVertices() != size_t(nLocalRows)) 00293 fail = 1; 00294 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalNumVertices", 1) 00295 00296 if (model->getGlobalNumVertices() != size_t(nGlobalRows)) 00297 fail = 1; 00298 TEST_FAIL_AND_EXIT(*comm, !fail, "getGlobalNumVertices", 1) 00299 00300 size_t num = (removeSelfEdges ? (nLocalNZ-numLocalDiags) : nLocalNZ); 00301 00302 if (model->getLocalNumGlobalEdges() != num) 00303 fail = 1; 00304 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalNumGlobalEdges", 1) 00305 00306 num = (removeSelfEdges ? (totalLocalNbors - numLocalDiags) : totalLocalNbors); 00307 00308 if (model->getLocalNumLocalEdges() != num) 00309 fail = 1; 00310 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalNumLocalEdges", 1) 00311 00312 num = (removeSelfEdges ? (nGlobalNZ-numGlobalDiags) : nGlobalNZ); 00313 00314 if (model->getGlobalNumEdges() != num) 00315 fail = 1; 00316 TEST_FAIL_AND_EXIT(*comm, !fail, "getGlobalNumEdges", 1) 00317 00318 if (model->getNumWeightsPerVertex() != nVtxWeights) 00319 fail = 1; 00320 TEST_FAIL_AND_EXIT(*comm, !fail, "getNumWeightsPerVertex", 1) 00321 00322 if (model->getNumWeightsPerEdge() != 0) 00323 fail = 1; 00324 TEST_FAIL_AND_EXIT(*comm, !fail, "getNumWeightsPerEdge", 1) 00325 00326 if (model->getCoordinateDim() != coordDim) 00327 fail = 1; 00328 TEST_FAIL_AND_EXIT(*comm, !fail, "getCoordinateDim", 1) 00329 00330 if (rank == 0) std::cout << " Checking vertices" << std::endl; 00331 ArrayView<const zgno_t> vertexGids; 00332 ArrayView<input_t> crds; 00333 ArrayView<input_t> wgts; 00334 00335 try{ 00336 model->getVertexList(vertexGids, crds, wgts); 00337 } 00338 catch (std::exception &e){ 00339 std::cerr << rank << ") Error " << e.what() << std::endl; 00340 fail = 1; 00341 } 00342 TEST_FAIL_AND_EXIT(*comm, !fail, "getVertexList", 1) 00343 00344 if (vertexGids.size() != nLocalRows) 00345 fail = 1; 00346 TEST_FAIL_AND_EXIT(*comm, !fail, "getVertexList size", 1) 00347 00348 // We know model stores things in same order we gave it. 00349 if (idsAreConsecutive){ 00350 zgno_t minLocalGID = rowMap->getMinGlobalIndex(); 00351 for (zlno_t i=0; i < nLocalRows; i++){ 00352 if (vertexGids[i] != minLocalGID + i) { 00353 fail = 1; 00354 break; 00355 } 00356 } 00357 } 00358 else{ // round robin ids 00359 if (consecutiveIdsRequested) { 00360 zgno_t myFirstRow; 00361 zgno_t tnLocalRows = nLocalRows; 00362 scan(*comm, Teuchos::REDUCE_SUM, 1, &tnLocalRows, &myFirstRow); 00363 myFirstRow -= nLocalRows; 00364 for (zlno_t i=0; i < nLocalRows; i++){ 00365 if (vertexGids[i] != myFirstRow+i){ 00366 std::cout << rank << " Row " << i << " of " << nLocalRows 00367 << " myFirstRow+i " << myFirstRow+i 00368 << " vertexGids " << vertexGids[i] 00369 << std::endl; 00370 fail = 1; 00371 break; 00372 } 00373 } 00374 } 00375 else { 00376 zgno_t myGid = rank; 00377 for (zlno_t i=0; i < nLocalRows; i++, myGid += nprocs){ 00378 if (vertexGids[i] != myGid){ 00379 std::cout << rank << " Row " << i << " of " << nLocalRows 00380 << " myGid " << myGid << " vertexGids " << vertexGids[i] 00381 << std::endl; 00382 fail = 1; 00383 break; 00384 } 00385 } 00386 } 00387 } 00388 00389 TEST_FAIL_AND_EXIT(*comm, !fail, "vertex gids", 1) 00390 00391 if ((crds.size() != coordDim) || (wgts.size() != nVtxWeights)) 00392 fail = 1; 00393 TEST_FAIL_AND_EXIT(*comm, !fail, "coord or weight array size", 1) 00394 00395 for (int i=0; !fail && i < coordDim; i++){ 00396 for (zlno_t j=0; j < nLocalRows; j++){ 00397 if (crds[i][j] != 100000*i + j){ 00398 fail = 1; 00399 break; 00400 } 00401 } 00402 } 00403 TEST_FAIL_AND_EXIT(*comm, !fail, "coord values", 1) 00404 00405 for (int i=0; !fail && i < nVtxWeights; i++){ 00406 if (nnzWgtIdx == i){ 00407 for (zlno_t j=0; j < nLocalRows; j++){ 00408 zscalar_t val = numNbors[j]; 00409 if (removeSelfEdges && haveDiag[j]) 00410 val -= 1; 00411 if (wgts[i][j] != val){ 00412 fail = 1; 00413 break; 00414 } 00415 } 00416 } 00417 else{ 00418 for (zlno_t j=0; j < nLocalRows; j++){ 00419 if (wgts[i][j] != 200000*i + j){ 00420 fail = 1; 00421 break; 00422 } 00423 } 00424 } 00425 } 00426 00427 TEST_FAIL_AND_EXIT(*comm, !fail, "row weight values", 1) 00428 00429 if (rank == 0) std::cout << " Checking edges" << std::endl; 00430 ArrayView<const zgno_t> edgeGids; 00431 ArrayView<const int> procIds; 00432 ArrayView<const zlno_t> offsets; 00433 size_t numEdges=0; 00434 00435 try{ 00436 numEdges = model->getEdgeList(edgeGids, procIds, offsets, wgts); 00437 } 00438 catch(std::exception &e){ 00439 std::cerr << rank << ") Error " << e.what() << std::endl; 00440 fail = 1; 00441 } 00442 TEST_FAIL_AND_EXIT(*comm, !fail, "getEdgeList", 1) 00443 00444 TEST_FAIL_AND_EXIT(*comm, wgts.size() == 0, "edge weights present", 1) 00445 00446 num = 0; 00447 for (ArrayView<const zlno_t>::size_type i=0; i < offsets.size()-1; i++){ 00448 size_t edgeListSize = offsets[i+1] - offsets[i]; 00449 num += edgeListSize; 00450 size_t val = numNbors[i]; 00451 if (removeSelfEdges && haveDiag[i]) 00452 val--; 00453 if (edgeListSize != val){ 00454 fail = 1; 00455 break; 00456 } 00457 } 00458 00459 TEST_FAIL_AND_EXIT(*comm, numEdges==num, "getEdgeList size", 1) 00460 00461 if (nGlobalRows < 200){ 00462 if (rank == 0) 00463 std::cout << "Printing global graph now " << nGlobalRows << std::endl; 00464 printGraph(nLocalRows, vertexGids.getRawPtr(), NULL, 00465 edgeGids.getRawPtr(), procIds.getRawPtr(), offsets.getRawPtr(), comm); 00466 } 00467 else{ 00468 if (rank==0) 00469 std::cout << " " << nGlobalRows << " total rows" << std::endl; 00470 } 00471 00472 // Get graph restricted to this process 00473 00474 if (rank == 0) std::cout << " Checking local edges" << std::endl; 00475 ArrayView<const zlno_t> localEdges; 00476 ArrayView<const zlno_t> localOffsets; 00477 size_t numLocalNeighbors=0; 00478 00479 try{ 00480 numLocalNeighbors= model->getLocalEdgeList(localEdges, localOffsets, wgts); 00481 } 00482 catch(std::exception &e){ 00483 std::cout << rank << ") Error " << e.what() << std::endl; 00484 std::cerr << rank << ") Error " << e.what() << std::endl; 00485 fail = 1; 00486 } 00487 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalEdgeList", 1) 00488 TEST_FAIL_AND_EXIT(*comm, ((localOffsets.size()-1) == nLocalRows), 00489 "getLocalEdgeList localOffsets.size", 1) 00490 00491 num = 0; 00492 for (zlno_t i=0; i < nLocalRows; i++){ 00493 size_t edgeListSize = localOffsets[i+1] - localOffsets[i]; 00494 num += edgeListSize; 00495 size_t val = numLocalNbors[i]; 00496 if (removeSelfEdges && haveDiag[i]) 00497 val--; 00498 if (edgeListSize != val){ 00499 std::cout << rank << "vtx " << i << " of " << localOffsets.size() 00500 << " Number of local edges in model " << edgeListSize 00501 << " not equal to expected number of edges " << val 00502 << std::endl; 00503 fail = 1; 00504 break; 00505 } 00506 } 00507 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalEdgeList list sizes", 1) 00508 00509 TEST_FAIL_AND_EXIT(*comm, numLocalNeighbors==num, 00510 "getLocalEdgeList sum size", 1) 00511 00512 fail = ((removeSelfEdges ? size_t(totalLocalNbors-numLocalDiags) 00513 : size_t(totalLocalNbors)) 00514 != numLocalNeighbors); 00515 TEST_FAIL_AND_EXIT(*comm, !fail, "getLocalEdgeList total size", 1) 00516 00517 if (nGlobalRows < 200){ 00518 if (totalLocalNbors == 0){ 00519 if (rank == 0) 00520 std::cout << " Graph of local edges is empty" << std::endl; 00521 } 00522 else{ 00523 printGraph(nLocalRows, vertexGids.getRawPtr(), 00524 localEdges.getRawPtr(), NULL, NULL, localOffsets.getRawPtr(), comm); 00525 } 00526 } 00527 00528 if (rank == 0) std::cout << " Cleaning up" << std::endl; 00529 delete model; 00530 00531 if (nLocalRows){ 00532 delete [] numNbors; 00533 delete [] numLocalNbors; 00534 delete [] haveDiag; 00535 00536 if (nVtxWeights > 0){ 00537 for (int i=0; i < nVtxWeights; i++){ 00538 if (rowWeights[i]) 00539 delete [] rowWeights[i]; 00540 } 00541 delete [] rowWeights; 00542 } 00543 00544 if (coordDim > 0){ 00545 delete via; 00546 delete [] gids; 00547 for (int i=0; i < coordDim; i++){ 00548 if (coords[i]) 00549 delete [] coords[i]; 00550 } 00551 delete [] coords; 00552 } 00553 } 00554 00555 if (rank==0) std::cout << " OK" << std::endl; 00556 } 00557 00559 void testGraphModel(string fname, zgno_t xdim, zgno_t ydim, zgno_t zdim, 00560 const RCP<const Comm<int> > &comm, 00561 int nVtxWeights, int nnzWgtIdx, int coordDim, 00562 bool consecutiveIdsRequested, bool removeSelfEdges) 00563 { 00564 int rank = comm->getRank(); 00565 00566 if (rank==0){ 00567 cout << endl << "=======================" << endl; 00568 if (fname.size() > 0) 00569 cout << endl << "Test parameters: file name " << fname << endl; 00570 else{ 00571 cout << endl << "Test parameters: dimension "; 00572 cout << xdim << "x" << ydim << "x" << zdim << endl; 00573 } 00574 00575 cout << "Num Vertex Weights: " << nVtxWeights << endl; 00576 if (nnzWgtIdx >= 0) 00577 cout << " Dimension " << nnzWgtIdx << " is number of neighbors" << endl; 00578 00579 cout << "Coordinate dim: " << coordDim << endl; 00580 cout << "Request consecutive vertex gids: "; 00581 cout << (consecutiveIdsRequested ? "yes" : "no") << endl; 00582 cout << "Request to remove self edges: "; 00583 cout << (removeSelfEdges ? "yes" : "no") << endl; 00584 } 00585 00586 typedef Tpetra::CrsMatrix<zscalar_t, zlno_t, zgno_t, znode_t> tcrsMatrix_t; 00587 00588 // Input generator 00589 UserInputForTests *uinput; 00590 00591 if (fname.size() > 0) 00592 uinput = new UserInputForTests(testDataFilePath, fname, comm, true); 00593 else 00594 uinput = new UserInputForTests(xdim,ydim,zdim,string(""), comm, true, true); 00595 00596 RCP<tcrsMatrix_t> M = uinput->getUITpetraCrsMatrix(); 00597 00598 // Row Ids of test input are already consecutive 00599 00600 RCP<const tcrsMatrix_t> Mconsec = rcp_const_cast<const tcrsMatrix_t>(M); 00601 00602 RCP<const Tpetra::CrsGraph<zlno_t, zgno_t> > graph = Mconsec->getCrsGraph(); 00603 00604 printTpetraGraph<zlno_t, zgno_t>(comm, *graph, cout, 100, 00605 "Graph with consecutive IDs"); 00606 00607 if (rank == 0) 00608 std::cout << " TEST MatrixAdapter with Consecutive IDs" << std::endl; 00609 bool idsAreConsecutive = true; 00610 00611 testAdapter<baseMAdapter_t,xMAdapter_t,tcrsMatrix_t>(Mconsec, graph, comm, 00612 idsAreConsecutive, 00613 nVtxWeights, 0, 00614 nnzWgtIdx, coordDim, 00615 consecutiveIdsRequested, 00616 removeSelfEdges); 00617 if (rank == 0) 00618 std::cout << " TEST GraphAdapter with Consecutive IDs" << std::endl; 00619 testAdapter<baseGAdapter_t,xGAdapter_t,tcrsGraph_t>(graph, graph, comm, 00620 idsAreConsecutive, 00621 nVtxWeights, 1, 00622 nnzWgtIdx, coordDim, 00623 consecutiveIdsRequested, 00624 removeSelfEdges); 00625 00626 // Do a round robin migration so that global IDs are not consecutive. 00627 00628 Array<zgno_t> myNewRows; 00629 int nprocs = comm->getSize(); 00630 for (size_t i=rank; i < Mconsec->getGlobalNumRows(); i+=nprocs) 00631 myNewRows.push_back(i); 00632 00633 RCP<const tcrsMatrix_t> Mnonconsec = 00634 Zoltan2::XpetraTraits<tcrsMatrix_t>::doMigration( 00635 Mconsec, myNewRows.size(), myNewRows.getRawPtr()); 00636 00637 graph = Mnonconsec->getCrsGraph(); 00638 00639 printTpetraGraph<zlno_t, zgno_t>(comm, *graph, cout, 100, 00640 "Graph with non-consecutive IDs"); 00641 00642 if (rank == 0) 00643 std::cout << " TEST MatrixAdapter with Round-Robin IDs" << std::endl; 00644 idsAreConsecutive = false; 00645 00646 testAdapter<baseMAdapter_t,xMAdapter_t,tcrsMatrix_t>(Mnonconsec, graph, comm, 00647 idsAreConsecutive, 00648 nVtxWeights, 0, 00649 nnzWgtIdx, coordDim, 00650 consecutiveIdsRequested, 00651 removeSelfEdges); 00652 00653 if (rank == 0) 00654 std::cout << " TEST GraphAdapter with Round-Robin IDs" << std::endl; 00655 testAdapter<baseGAdapter_t,xGAdapter_t,tcrsGraph_t>(graph, graph, comm, 00656 idsAreConsecutive, 00657 nVtxWeights, 0, 00658 nnzWgtIdx, coordDim, 00659 consecutiveIdsRequested, 00660 removeSelfEdges); 00661 00662 delete uinput; 00663 } 00664 00666 int main(int argc, char *argv[]) 00667 { 00668 Teuchos::GlobalMPISession session(&argc, &argv); 00669 Teuchos::RCP<const Teuchos::Comm<int> > comm = 00670 Teuchos::DefaultComm<int>::getComm(); 00671 00672 int rank = comm->getRank(); 00673 00674 int nVtxWeights=0; 00675 int nnzWgtIdx = -1; 00676 int coordDim=0; 00677 bool consecutiveIdsRequested=false, removeSelfEdges=false; 00678 string fname("simple"); 00679 00680 if (rank == 0) 00681 std::cout << "TESTING base case" << std::endl; 00682 testGraphModel(fname, 0, 0, 0, comm, 00683 nVtxWeights, nnzWgtIdx, coordDim, 00684 consecutiveIdsRequested, removeSelfEdges); 00685 00686 if (rank == 0) 00687 std::cout << "TESTING with row weights" << std::endl; 00688 nVtxWeights = 1; 00689 00690 testGraphModel(fname, 0, 0, 0, comm, 00691 nVtxWeights, nnzWgtIdx, coordDim, 00692 consecutiveIdsRequested, removeSelfEdges); 00693 00694 if (rank == 0) 00695 std::cout << "TESTING with weights = nnz" << std::endl; 00696 nnzWgtIdx = 1; 00697 00698 testGraphModel(fname, 0, 0, 0, comm, 00699 nVtxWeights, nnzWgtIdx, coordDim, 00700 consecutiveIdsRequested, removeSelfEdges); 00701 00702 if (rank == 0) 00703 std::cout << "TESTING with multiple row weights and coords" << std::endl; 00704 nVtxWeights = 2; 00705 coordDim = 3; 00706 00707 testGraphModel(fname, 0, 0, 0, comm, 00708 nVtxWeights, nnzWgtIdx, coordDim, 00709 consecutiveIdsRequested, removeSelfEdges); 00710 00711 if (rank == 0) 00712 std::cout << "TESTING with consecutiveIdsRequested" << std::endl; 00713 consecutiveIdsRequested = true; 00714 00715 testGraphModel(fname, 0, 0, 0, comm, 00716 nVtxWeights, nnzWgtIdx, coordDim, 00717 consecutiveIdsRequested, removeSelfEdges); 00718 00719 if (rank == 0) 00720 std::cout << "TESTING with removeSelfEdges" << std::endl; 00721 removeSelfEdges = true; 00722 00723 testGraphModel(fname, 0, 0, 0, comm, 00724 nVtxWeights, nnzWgtIdx, coordDim, 00725 consecutiveIdsRequested, removeSelfEdges); 00726 00727 if (rank == 0) 00728 std::cout << "TESTING with consecutiveIdsRequested=false" << std::endl; 00729 consecutiveIdsRequested = false; 00730 00731 testGraphModel(fname, 0, 0, 0, comm, 00732 nVtxWeights, nnzWgtIdx, coordDim, 00733 consecutiveIdsRequested, removeSelfEdges); 00734 00735 if (rank==0) 00736 cout << "PASS" << endl; 00737 00738 return 0; 00739 } 00740
1.7.6.1