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