Zoltan2
XpetraVectorInput.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 
00052 #include <string>
00053 
00054 #include <Zoltan2_XpetraMultiVectorAdapter.hpp>
00055 #include <Zoltan2_InputTraits.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 
00071 typedef UserInputForTests uinput_t;
00072 typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> tvector_t;
00073 typedef Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> xvector_t;
00074 typedef Epetra_Vector evector_t;
00075 
00076 void printVector(RCP<const Comm<int> > &comm, zlno_t vlen,
00077     const zgno_t *vtxIds, const zscalar_t *vals)
00078 {
00079   int rank = comm->getRank();
00080   int nprocs = comm->getSize();
00081   comm->barrier();
00082   for (int p=0; p < nprocs; p++){
00083     if (p == rank){
00084       std::cout << rank << ":" << std::endl;
00085       for (zlno_t i=0; i < vlen; i++){
00086         std::cout << " " << vtxIds[i] << ": " << vals[i] << std::endl;
00087       }
00088       std::cout.flush();
00089     }
00090     comm->barrier();
00091   }
00092   comm->barrier();
00093 }
00094 
00095 template <typename User>
00096 int verifyInputAdapter(
00097   Zoltan2::XpetraMultiVectorAdapter<User> &ia, tvector_t &vector, int wdim, 
00098     zscalar_t **weights, int *strides)
00099 {
00100   RCP<const Comm<int> > comm = vector.getMap()->getComm();
00101   int fail = 0, gfail=0;
00102 
00103   if (!fail && ia.getNumEntriesPerID() !=1) 
00104     fail = 42;
00105 
00106   if (!fail && ia.getNumWeightsPerID() !=wdim) 
00107     fail = 41;
00108 
00109   if (!fail && ia.getLocalNumIDs() != vector.getLocalLength())
00110     fail = 4;
00111 
00112   gfail = globalFail(comm, fail);
00113 
00114   if (!gfail){
00115     const zgno_t *vtxIds=NULL;
00116     const zscalar_t *vals=NULL;
00117     int stride;
00118 
00119     size_t nvals = ia.getLocalNumIDs();
00120     if (nvals != vector.getLocalLength())
00121       fail = 8;
00122 
00123     ia.getIDsView(vtxIds);
00124     ia.getEntriesView(vals, stride);
00125 
00126     if (!fail && stride != 1)
00127       fail = 10;
00128 
00129     gfail = globalFail(comm, fail);
00130 
00131     if (gfail == 0){
00132       printVector(comm, nvals, vtxIds, vals);
00133     }
00134   }
00135 
00136   if (!gfail && wdim){
00137     const zscalar_t *wgt =NULL;
00138     int stride;
00139 
00140     for (int w=0; !fail && w < wdim; w++){
00141       ia.getWeightsView(wgt, stride, w);
00142 
00143       if (!fail && stride != strides[w])
00144         fail = 101;
00145 
00146       for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
00147         if (wgt[v*stride] != weights[w][v*stride])
00148           fail=102;
00149       }
00150     }
00151 
00152     gfail = globalFail(comm, fail);
00153   }
00154 
00155   return gfail;
00156 }
00157 
00158 int main(int argc, char *argv[])
00159 {
00160   Teuchos::GlobalMPISession session(&argc, &argv);
00161   RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
00162   int rank = comm->getRank();
00163   int fail = 0, gfail=0;
00164 
00165   // Create object that can give us test Tpetra, Xpetra
00166   // and Epetra vectors for testing.
00167 
00168   RCP<uinput_t> uinput;
00169 
00170   try{
00171     uinput = 
00172       rcp(new uinput_t(testDataFilePath,std::string("simple"), comm, true));
00173   }
00174   catch(std::exception &e){
00175     TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
00176   }
00177 
00178   RCP<tvector_t> tV;     // original vector (for checking)
00179   RCP<tvector_t> newV;   // migrated vector
00180 
00181   tV = uinput->getUITpetraVector();
00182   size_t vlen = tV->getLocalLength();
00183   Teuchos::ArrayView<const zgno_t> rowGids = tV->getMap()->getNodeElementList();
00184 
00185   // To test migration in the input adapter we need a Solution
00186   // object.  The Solution needs an IdentifierMap.
00187 
00188   typedef Zoltan2::IdentifierMap<tvector_t> idmap_t;
00189 
00190   RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment);
00191 
00192   ArrayRCP<const zgno_t> gidArray = arcpFromArrayView(rowGids);
00193   RCP<const idmap_t> idMap = rcp(new idmap_t(env, comm, gidArray));
00194 
00195   int nWeights = 1;
00196 
00197   typedef Zoltan2::XpetraMultiVectorAdapter<tvector_t> adapter_t;
00198   typedef adapter_t::part_t part_t;
00199 
00200   part_t *p = new part_t [vlen];
00201   memset(p, 0, sizeof(part_t) * vlen);
00202   ArrayRCP<part_t> solnParts(p, 0, vlen, true);
00203 
00204   std::vector<const zscalar_t *> emptyWeights;
00205   std::vector<int> emptyStrides;
00206 
00207   Zoltan2::PartitioningSolution<adapter_t> solution(
00208     env, comm, idMap, nWeights);
00209   solution.setParts(gidArray, solnParts, true);
00210 
00212   // User object is Tpetra::Vector, no weights
00213   if (!gfail){ 
00214     RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
00215     RCP<adapter_t> tVInput;
00216   
00217     try {
00218       tVInput = 
00219         rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(ctV, 
00220           emptyWeights, emptyStrides));
00221     }
00222     catch (std::exception &e){
00223       TEST_FAIL_AND_EXIT(*comm, 0, 
00224         string("XpetraMultiVectorAdapter ")+e.what(), 1);
00225     }
00226   
00227     if (rank==0){
00228       std::cout << "Constructed with ";
00229       std::cout  << "Tpetra::Vector" << std::endl;
00230     }
00231     
00232     fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
00233   
00234     gfail = globalFail(comm, fail);
00235   
00236     if (!gfail){
00237       tvector_t *vMigrate = NULL;
00238       try{
00239         tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
00240         newV = rcp(vMigrate);
00241       }
00242       catch (std::exception &e){
00243         fail = 11;
00244       }
00245 
00246       gfail = globalFail(comm, fail);
00247   
00248       if (!gfail){
00249         RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
00250         RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
00251         try{
00252           newInput = rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(cnewV,
00253             emptyWeights, emptyStrides));
00254         }
00255         catch (std::exception &e){
00256           TEST_FAIL_AND_EXIT(*comm, 0, 
00257             string("XpetraMultiVectorAdapter 2 ")+e.what(), 1);
00258         }
00259   
00260         if (rank==0){
00261           std::cout << "Constructed with ";
00262           std::cout << "Tpetra::Vector migrated to proc 0" << std::endl;
00263         }
00264         fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
00265         if (fail) fail += 100;
00266         gfail = globalFail(comm, fail);
00267       }
00268     }
00269     if (gfail){
00270       printFailureCode(comm, fail);
00271     }
00272   }
00273 
00275   // User object is Xpetra::Vector
00276   if (!gfail){ 
00277     RCP<xvector_t> xV = uinput->getUIXpetraVector();
00278     RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
00279     RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
00280   
00281     try {
00282       xVInput = 
00283         rcp(new Zoltan2::XpetraMultiVectorAdapter<xvector_t>(cxV,
00284           emptyWeights, emptyStrides));
00285     }
00286     catch (std::exception &e){
00287       TEST_FAIL_AND_EXIT(*comm, 0, 
00288         string("XpetraMultiVectorAdapter 3 ")+e.what(), 1);
00289     }
00290   
00291     if (rank==0){
00292       std::cout << "Constructed with ";
00293       std::cout << "Xpetra::Vector" << std::endl;
00294     }
00295     fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
00296   
00297     gfail = globalFail(comm, fail);
00298   
00299     if (!gfail){
00300       xvector_t *vMigrate =NULL;
00301       try{
00302         xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
00303       }
00304       catch (std::exception &e){
00305         fail = 11;
00306       }
00307   
00308       gfail = globalFail(comm, fail);
00309   
00310       if (!gfail){
00311         RCP<const xvector_t> cnewV(vMigrate);
00312         RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
00313         try{
00314           newInput = 
00315             rcp(new Zoltan2::XpetraMultiVectorAdapter<xvector_t>(cnewV, 
00316               emptyWeights, emptyStrides));
00317         }
00318         catch (std::exception &e){
00319           TEST_FAIL_AND_EXIT(*comm, 0, 
00320             string("XpetraMultiVectorAdapter 4 ")+e.what(), 1);
00321         }
00322   
00323         if (rank==0){
00324           std::cout << "Constructed with ";
00325           std::cout << "Xpetra::Vector migrated to proc 0" << std::endl;
00326         }
00327         fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
00328         if (fail) fail += 100;
00329         gfail = globalFail(comm, fail);
00330       }
00331     }
00332     if (gfail){
00333       printFailureCode(comm, fail);
00334     }
00335   }
00336 
00337 #ifdef HAVE_EPETRA_DATA_TYPES
00338 
00339   // User object is Epetra_Vector
00340   if (!gfail){ 
00341     RCP<evector_t> eV = uinput->getUIEpetraVector();
00342     RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
00343     RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
00344   
00345     try {
00346       eVInput = 
00347         rcp(new Zoltan2::XpetraMultiVectorAdapter<evector_t>(ceV,
00348           emptyWeights, emptyStrides));
00349     }
00350     catch (std::exception &e){
00351       TEST_FAIL_AND_EXIT(*comm, 0, 
00352         string("XpetraMultiVectorAdapter 5 ")+e.what(), 1);
00353     }
00354   
00355     if (rank==0){
00356       std::cout << "Constructed with ";
00357       std::cout << "Epetra_Vector" << std::endl;
00358     }
00359     fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
00360   
00361     gfail = globalFail(comm, fail);
00362   
00363     if (!gfail){
00364       evector_t *vMigrate =NULL;
00365       try{
00366         eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
00367       }
00368       catch (std::exception &e){
00369         fail = 11;
00370       }
00371   
00372       gfail = globalFail(comm, fail);
00373   
00374       if (!gfail){
00375         RCP<const evector_t> cnewV(vMigrate, true);
00376         RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
00377         try{
00378           newInput = 
00379             rcp(new Zoltan2::XpetraMultiVectorAdapter<evector_t>(cnewV, 
00380               emptyWeights, emptyStrides));
00381         }
00382         catch (std::exception &e){
00383           TEST_FAIL_AND_EXIT(*comm, 0, 
00384             string("XpetraMultiVectorAdapter 6 ")+e.what(), 1);
00385         }
00386   
00387         if (rank==0){
00388            std::cout << "Constructed with ";
00389            std::cout << "Epetra_Vector migrated to proc 0" << std::endl;
00390         }
00391         fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
00392         if (fail) fail += 100;
00393         gfail = globalFail(comm, fail);
00394       }
00395     }
00396     if (gfail){
00397       printFailureCode(comm, fail);
00398     }
00399   }
00400 #endif
00401 
00403   // DONE
00404 
00405   if (rank==0)
00406     std::cout << "PASS" << std::endl;
00407 }