Zoltan2
XpetraCrsGraphInput.cpp
Go to the documentation of this file.
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 // Basic testing of Zoltan2::XpetraCrsGraphAdapter
00052 #include <string>
00053 
00054 #include <Zoltan2_XpetraCrsGraphAdapter.hpp>
00055 #include <Zoltan2_PartitioningSolution.hpp>
00056 #include <Zoltan2_TestHelpers.hpp>
00057 
00058 #include <Teuchos_GlobalMPISession.hpp>
00059 #include <Teuchos_DefaultComm.hpp>
00060 #include <Teuchos_RCP.hpp>
00061 #include <Teuchos_Comm.hpp>
00062 #include <Teuchos_CommHelpers.hpp>
00063 
00064 using namespace std;
00065 using Teuchos::RCP;
00066 using Teuchos::rcp;
00067 using Teuchos::rcp_const_cast;
00068 using Teuchos::Comm;
00069 using Teuchos::DefaultComm;
00070 using Teuchos::Array;
00071 using Teuchos::ArrayView;
00072 
00073 typedef UserInputForTests uinput_t;
00074 typedef Tpetra::CrsGraph<zlno_t, zgno_t, znode_t> tgraph_t;
00075 typedef Xpetra::CrsGraph<zlno_t, zgno_t, znode_t> xgraph_t;
00076 typedef Epetra_CrsGraph egraph_t;
00077 
00078 void printGraph(RCP<const Comm<int> > &comm, zlno_t nvtx,
00079     const zgno_t *vtxIds, const zlno_t *offsets, const zgno_t *edgeIds)
00080 {
00081   int rank = comm->getRank();
00082   int nprocs = comm->getSize();
00083   comm->barrier();
00084   for (int p=0; p < nprocs; p++){
00085     if (p == rank){
00086       std::cout << rank << ":" << std::endl;
00087       for (zlno_t i=0; i < nvtx; i++){
00088         std::cout << " vertex " << vtxIds[i] << ": ";
00089         for (zlno_t j=offsets[i]; j < offsets[i+1]; j++){
00090           std::cout << edgeIds[j] << " ";
00091         }
00092         std::cout << std::endl;
00093       }
00094       std::cout.flush();
00095     }
00096     comm->barrier();
00097   }
00098   comm->barrier();
00099 }
00100 
00101 template <typename User>
00102 int verifyInputAdapter(
00103   Zoltan2::XpetraCrsGraphAdapter<User> &ia,
00104   tgraph_t &graph
00105 )
00106 {
00107   RCP<const Comm<int> > comm = graph.getComm();
00108   int fail = 0, gfail=0;
00109 
00110   if (!fail && 
00111       ia.getLocalNumVertices() != graph.getNodeNumRows())
00112     fail = 4;
00113 
00114   if (!fail && 
00115       ia.getLocalNumEdges() != graph.getNodeNumEntries())
00116       fail = 6;
00117 
00118   gfail = globalFail(comm, fail);
00119 
00120   const zgno_t *vtxIds=NULL, *edgeIds=NULL;
00121   const zlno_t *offsets=NULL;
00122   size_t nvtx=0;
00123 
00124   if (!gfail){
00125 
00126     nvtx = ia.getLocalNumVertices();
00127     ia.getVertexIDsView(vtxIds);
00128     ia.getEdgesView(offsets, edgeIds);
00129 
00130     if (nvtx != graph.getNodeNumRows())
00131       fail = 8;
00132 
00133     gfail = globalFail(comm, fail);
00134 
00135     if (gfail == 0){
00136       printGraph(comm, nvtx, vtxIds, offsets, edgeIds);
00137     }
00138     else{
00139       if (!fail) fail = 10;
00140     }
00141   }
00142   return fail;
00143 }
00144 
00145 int main(int argc, char *argv[])
00146 {
00147   Teuchos::GlobalMPISession session(&argc, &argv);
00148   RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
00149   int rank = comm->getRank();
00150   int fail = 0, gfail=0;
00151 
00152   // Create an object that can give us test Tpetra, Xpetra
00153   // and Epetra graphs for testing.
00154 
00155   RCP<uinput_t> uinput;
00156 
00157   try{
00158     uinput =
00159       rcp(new uinput_t(testDataFilePath,std::string("simple"), comm, true));
00160   }
00161   catch(std::exception &e){
00162     TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
00163   }
00164 
00165   RCP<tgraph_t> tG;     // original graph (for checking)
00166   RCP<tgraph_t> newG;   // migrated graph
00167 
00168   tG = uinput->getUITpetraCrsGraph();
00169   size_t nvtx = tG->getNodeNumRows();
00170   ArrayView<const zgno_t> rowGids = tG->getRowMap()->getNodeElementList();
00171 
00172   // To test migration in the input adapter we need a Solution
00173   // object.  The Solution needs an IdentifierMap.
00174   // Our solution just assigns all objects to part zero.
00175 
00176   typedef Zoltan2::IdentifierMap<tgraph_t> idmap_t;
00177 
00178   RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment);
00179 
00180   ArrayRCP<const zgno_t> gidArray = arcpFromArrayView(rowGids);
00181   RCP<const idmap_t> idMap = rcp(new idmap_t(env, comm, gidArray));
00182 
00183   int nWeights = 1;
00184 
00185   typedef Zoltan2::XpetraCrsGraphAdapter<tgraph_t>  adapter_t;
00186   typedef Zoltan2::PartitioningSolution<adapter_t> soln_t;
00187   typedef adapter_t::part_t part_t;
00188 
00189   part_t *p = new part_t [nvtx];
00190   memset(p, 0, sizeof(part_t) * nvtx);
00191   ArrayRCP<part_t> solnParts(p, 0, nvtx, true);
00192 
00193   soln_t solution(env, comm, idMap, nWeights);
00194   solution.setParts(gidArray, solnParts, true);
00195 
00197   // User object is Tpetra::CrsGraph
00198   if (!gfail){
00199     RCP<const tgraph_t> ctG = rcp_const_cast<const tgraph_t>(tG);
00200     RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > tGInput;
00201 
00202     try {
00203       tGInput =
00204         rcp(new Zoltan2::XpetraCrsGraphAdapter<tgraph_t>(ctG));
00205     }
00206     catch (std::exception &e){
00207       TEST_FAIL_AND_EXIT(*comm, 0,
00208         string("XpetraCrsGraphAdapter ")+e.what(), 1);
00209     }
00210 
00211     if (rank==0)
00212       std::cout << "Input adapter for Tpetra::CrsGraph" << std::endl;
00213 
00214     fail = verifyInputAdapter<tgraph_t>(*tGInput, *tG);
00215 
00216     gfail = globalFail(comm, fail);
00217 
00218     if (!gfail){
00219       tgraph_t *mMigrate = NULL;
00220       try{
00221         tGInput->applyPartitioningSolution( *tG, mMigrate, solution);
00222         newG = rcp(mMigrate);
00223       }
00224       catch (std::exception &e){
00225         fail = 11;
00226       }
00227 
00228       gfail = globalFail(comm, fail);
00229 
00230       if (!gfail){
00231         RCP<const tgraph_t> cnewG = rcp_const_cast<const tgraph_t>(newG);
00232         RCP<Zoltan2::XpetraCrsGraphAdapter<tgraph_t> > newInput;
00233         try{
00234           newInput = rcp(new Zoltan2::XpetraCrsGraphAdapter<tgraph_t>(cnewG));
00235         }
00236         catch (std::exception &e){
00237           TEST_FAIL_AND_EXIT(*comm, 0,
00238             string("XpetraCrsGraphAdapter 2 ")+e.what(), 1);
00239         }
00240 
00241         if (rank==0){
00242           std::cout <<
00243            "Input adapter for Tpetra::CrsGraph migrated to proc 0" <<
00244            std::endl;
00245         }
00246         fail = verifyInputAdapter<tgraph_t>(*newInput, *newG);
00247         if (fail) fail += 100;
00248         gfail = globalFail(comm, fail);
00249       }
00250     }
00251     if (gfail){
00252       printFailureCode(comm, fail);
00253     }
00254   }
00255 
00257   // User object is Xpetra::CrsGraph
00258   if (!gfail){
00259     RCP<xgraph_t> xG = uinput->getUIXpetraCrsGraph();
00260     RCP<const xgraph_t> cxG = rcp_const_cast<const xgraph_t>(xG);
00261     RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > xGInput;
00262 
00263     try {
00264       xGInput =
00265         rcp(new Zoltan2::XpetraCrsGraphAdapter<xgraph_t>(cxG));
00266     }
00267     catch (std::exception &e){
00268       TEST_FAIL_AND_EXIT(*comm, 0,
00269         string("XpetraCrsGraphAdapter 3 ")+e.what(), 1);
00270     }
00271 
00272     if (rank==0){
00273       std::cout << "Input adapter for Xpetra::CrsGraph" << std::endl;
00274     }
00275     fail = verifyInputAdapter<xgraph_t>(*xGInput, *tG);
00276 
00277     gfail = globalFail(comm, fail);
00278 
00279     if (!gfail){
00280       xgraph_t *mMigrate =NULL;
00281       try{
00282         xGInput->applyPartitioningSolution( *xG, mMigrate, solution);
00283       }
00284       catch (std::exception &e){
00285         fail = 11;
00286       }
00287 
00288       gfail = globalFail(comm, fail);
00289 
00290       if (!gfail){
00291         RCP<const xgraph_t> cnewG(mMigrate);
00292         RCP<Zoltan2::XpetraCrsGraphAdapter<xgraph_t> > newInput;
00293         try{
00294           newInput =
00295             rcp(new Zoltan2::XpetraCrsGraphAdapter<xgraph_t>(cnewG));
00296         }
00297         catch (std::exception &e){
00298           TEST_FAIL_AND_EXIT(*comm, 0,
00299             string("XpetraCrsGraphAdapter 4 ")+e.what(), 1);
00300         }
00301 
00302         if (rank==0){
00303           std::cout <<
00304            "Input adapter for Xpetra::CrsGraph migrated to proc 0" <<
00305            std::endl;
00306         }
00307         fail = verifyInputAdapter<xgraph_t>(*newInput, *newG);
00308         if (fail) fail += 100;
00309         gfail = globalFail(comm, fail);
00310       }
00311     }
00312     if (gfail){
00313       printFailureCode(comm, fail);
00314     }
00315   }
00316 
00317 #ifdef HAVE_EPETRA_DATA_TYPES
00318 
00319   // User object is Epetra_CrsGraph
00320   if (!gfail){
00321     RCP<egraph_t> eG = uinput->getUIEpetraCrsGraph();
00322     RCP<const egraph_t> ceG = rcp_const_cast<const egraph_t>(eG);
00323     RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > eGInput;
00324 
00325     try {
00326       eGInput =
00327         rcp(new Zoltan2::XpetraCrsGraphAdapter<egraph_t>(ceG));
00328     }
00329     catch (std::exception &e){
00330       TEST_FAIL_AND_EXIT(*comm, 0,
00331         string("XpetraCrsGraphAdapter 5 ")+e.what(), 1);
00332     }
00333 
00334     if (rank==0){
00335       std::cout << "Input adapter for Epetra_CrsGraph" << std::endl;
00336     }
00337     fail = verifyInputAdapter<egraph_t>(*eGInput, *tG);
00338 
00339     gfail = globalFail(comm, fail);
00340 
00341     if (!gfail){
00342       egraph_t *mMigrate =NULL;
00343       try{
00344         eGInput->applyPartitioningSolution( *eG, mMigrate, solution);
00345       }
00346       catch (std::exception &e){
00347         fail = 11;
00348       }
00349 
00350       gfail = globalFail(comm, fail);
00351 
00352       if (!gfail){
00353         RCP<const egraph_t> cnewG(mMigrate, true);
00354         RCP<Zoltan2::XpetraCrsGraphAdapter<egraph_t> > newInput;
00355         try{
00356           newInput =
00357             rcp(new Zoltan2::XpetraCrsGraphAdapter<egraph_t>(cnewG));
00358         }
00359         catch (std::exception &e){
00360           TEST_FAIL_AND_EXIT(*comm, 0,
00361             string("XpetraCrsGraphAdapter 6 ")+e.what(), 1);
00362         }
00363 
00364         if (rank==0){
00365           std::cout <<
00366            "Input adapter for Epetra_CrsGraph migrated to proc 0" <<
00367            std::endl;
00368         }
00369         fail = verifyInputAdapter<egraph_t>(*newInput, *newG);
00370         if (fail) fail += 100;
00371         gfail = globalFail(comm, fail);
00372       }
00373     }
00374     if (gfail){
00375       printFailureCode(comm, fail);
00376     }
00377   }
00378 #endif
00379 
00381   // DONE
00382 
00383   if (rank==0)
00384     std::cout << "PASS" << std::endl;
00385 }