|
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 // Test for Zoltan2::BasicVectorAdapter 00047 00048 #include <Zoltan2_BasicVectorAdapter.hpp> 00049 #include <Zoltan2_TestHelpers.hpp> 00050 00051 #include <Teuchos_GlobalMPISession.hpp> 00052 #include <Teuchos_DefaultComm.hpp> 00053 #include <Teuchos_RCP.hpp> 00054 #include <Teuchos_CommHelpers.hpp> 00055 00056 using Teuchos::RCP; 00057 using Teuchos::Comm; 00058 using Teuchos::DefaultComm; 00059 00060 typedef Zoltan2::BasicUserTypes<zscalar_t, zzgid_t, zlno_t, zgno_t> userTypes_t; 00061 00062 int checkBasicVector( 00063 Zoltan2::BasicVectorAdapter<userTypes_t> *ia, int len, int glen, 00064 zzgid_t *ids, int mvdim, const zscalar_t **values, int *valueStrides, 00065 int wdim, const zscalar_t **weights, int *weightStrides) 00066 { 00067 int fail = 0; 00068 bool strideOne = false; 00069 00070 if (valueStrides == NULL) strideOne = true; 00071 00072 if (ia->getNumEntriesPerID() != mvdim) 00073 fail = 100; 00074 00075 if (!fail && ia->getNumWeightsPerID() != wdim) 00076 fail = 101; 00077 00078 if (!fail && ia->getLocalNumIDs() != size_t(len)) 00079 fail = 102; 00080 00081 const zzgid_t *idList; 00082 ia->getIDsView(idList); 00083 for (int i=0; !fail && i < len; i++) 00084 if (!fail && idList[i] != ids[i]) 00085 fail = 107; 00086 00087 for (int v=0; !fail && v < mvdim; v++){ 00088 const zscalar_t *vals; 00089 int correctStride = (strideOne ? 1 : valueStrides[v]); 00090 int stride; 00091 00092 ia->getEntriesView(vals, stride, v); 00093 00094 if (!fail && stride != correctStride) 00095 fail = 105; 00096 00097 for (int i=0; !fail && i < len; i++){ 00098 // TODO fix values check 00099 // if (vals[stride*i] != values[v][correctStride*i]) 00100 // fail = 106; 00101 } 00102 } 00103 00104 for (int w=0; !fail && w < wdim; w++){ 00105 const zscalar_t *wgts; 00106 int stride; 00107 00108 ia->getWeightsView(wgts, stride, w); 00109 00110 if (!fail && stride != weightStrides[w]) 00111 fail = 109; 00112 00113 for (int i=0; !fail && i < len; i++){ 00114 if (wgts[stride*i] != weights[w][weightStrides[w]*i]) 00115 fail = 110; 00116 } 00117 } 00118 00119 return fail; 00120 } 00121 00122 00123 int main(int argc, char *argv[]) 00124 { 00125 Teuchos::GlobalMPISession session(&argc, &argv); 00126 RCP<const Comm<int> > comm = DefaultComm<int>::getComm(); 00127 int rank = comm->getRank(); 00128 int nprocs = comm->getSize(); 00129 int fail = 0; 00130 00131 // Create a single vector and a strided multi-vector with 00132 // strided multi-weights. 00133 00134 zlno_t numLocalIds = 10; 00135 zzgid_t *myIds = new zzgid_t [numLocalIds]; 00136 zgno_t base = rank * numLocalIds; 00137 00138 int wdim = 2; 00139 zscalar_t *weights = new zscalar_t [numLocalIds*wdim]; 00140 int *weightStrides = new int [wdim]; 00141 const zscalar_t **weightPtrs = new const zscalar_t * [wdim]; 00142 00143 int mvdim = 3; 00144 zscalar_t *v_values= new zscalar_t [numLocalIds]; 00145 zscalar_t *mv_values= new zscalar_t [mvdim * numLocalIds]; 00146 int *valueStrides = new int [mvdim]; 00147 const zscalar_t **valuePtrs = new const zscalar_t * [mvdim]; 00148 00149 for (zlno_t i=0; i < numLocalIds; i++){ 00150 myIds[i] = base+i; 00151 00152 for (int w=0; w < wdim; w++) 00153 weights[w*numLocalIds + i] = w + 1 + nprocs - rank; 00154 00155 v_values[i] = numLocalIds-i; 00156 00157 for (int v=0; v < mvdim; v++) 00158 mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1); 00159 } 00160 00161 for (int w=0; w < wdim; w++){ 00162 weightStrides[w] = 1; 00163 weightPtrs[w] = weights + numLocalIds*w; 00164 } 00165 00166 for (int v=0; v < mvdim; v++){ 00167 valueStrides[v] = mvdim; 00168 valuePtrs[v] = mv_values + v; 00169 } 00170 00171 Zoltan2::BasicVectorAdapter<userTypes_t> *ia = NULL; 00172 00173 { 00174 // A Zoltan2::BasicVectorAdapter object with one vector and no weights 00175 00176 std::vector<const zscalar_t *> weightValues; 00177 std::vector<int> strides; 00178 00179 try{ 00180 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds, 00181 v_values, 1); 00182 } 00183 catch (std::exception &e){ 00184 fail = 1; 00185 } 00186 00187 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail); 00188 00189 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs, 00190 myIds, 1, valuePtrs, NULL, 0, NULL, NULL); 00191 00192 delete ia; 00193 00194 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 1", fail); 00195 } 00196 00197 { 00198 // A Zoltan2::BasicVectorAdapter object with one vector and weights 00199 00200 std::vector<const zscalar_t *> weightValues; 00201 std::vector<int> strides; 00202 00203 weightValues.push_back(weightPtrs[0]); 00204 weightValues.push_back(weightPtrs[1]); 00205 strides.push_back(1); 00206 strides.push_back(1); 00207 00208 try{ 00209 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds, 00210 v_values, 1, true, weightPtrs[0], 1); 00211 } 00212 catch (std::exception &e){ 00213 fail = 1; 00214 } 00215 00216 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail); 00217 00218 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs, 00219 myIds, 1, valuePtrs, NULL, 1, weightPtrs, weightStrides); 00220 00221 delete ia; 00222 00223 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 2", fail); 00224 } 00225 00226 { 00227 // A Zoltan2::BasicVectorAdapter object with a multivector and no weights 00228 00229 std::vector<const zscalar_t *> weightValues, values; 00230 std::vector<int> wstrides, vstrides; 00231 00232 for (int dim=0; dim < mvdim; dim++){ 00233 values.push_back(valuePtrs[dim]); 00234 vstrides.push_back(mvdim); 00235 } 00236 00237 00238 try{ 00239 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>( 00240 numLocalIds, myIds, values, vstrides, weightValues, wstrides); 00241 } 00242 catch (std::exception &e){ 00243 fail = 1; 00244 } 00245 00246 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail); 00247 00248 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs, 00249 myIds, mvdim, valuePtrs, valueStrides, 0, NULL, NULL); 00250 00251 delete ia; 00252 00253 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 3", fail); 00254 } 00255 00256 { 00257 // A Zoltan2::BasicVectorAdapter object with a multivector with weights 00258 00259 std::vector<const zscalar_t *> weightValues, values; 00260 std::vector<int> wstrides, vstrides; 00261 00262 for (int dim=0; dim < wdim; dim++){ 00263 weightValues.push_back(weightPtrs[dim]); 00264 wstrides.push_back(1); 00265 } 00266 00267 for (int dim=0; dim < mvdim; dim++){ 00268 values.push_back(valuePtrs[dim]); 00269 vstrides.push_back(mvdim); 00270 } 00271 00272 try{ 00273 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>( 00274 numLocalIds, myIds, values, vstrides, weightValues, wstrides); 00275 00276 } 00277 catch (std::exception &e){ 00278 fail = 1; 00279 } 00280 00281 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail); 00282 00283 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs, 00284 myIds, mvdim, valuePtrs, valueStrides, 00285 wdim, weightPtrs, weightStrides); 00286 00287 delete ia; 00288 00289 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 4", fail); 00290 } 00291 00292 if (rank == 0) 00293 std::cout << "PASS" << std::endl; 00294 00295 delete [] myIds; 00296 delete [] weights; 00297 delete [] weightStrides; 00298 delete [] weightPtrs; 00299 delete [] v_values; 00300 delete [] mv_values; 00301 delete [] valueStrides; 00302 delete [] valuePtrs; 00303 00304 return fail; 00305 } 00306
1.7.6.1