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