00001
00002 #include "mpi.h"
00003 #include <iostream>
00004 #include <cstdlib>
00005 #include <ctime>
00006 #include "petscvec.h"
00007 #include "dtypes.h"
00008 #include "omg.h"
00009 #include "parUtils.h"
00010 #include "externVars.h"
00011 #include "dendro.h"
00012
00013 int main(int argc, char ** argv ) {
00014
00015 int rank, npes;
00016
00017 PetscInitialize(&argc,&argv,NULL,NULL);
00018
00019 MPI_Comm_rank(MPI_COMM_WORLD,&rank);
00020 MPI_Comm_size(MPI_COMM_WORLD,&npes);
00021
00022
00023 {
00024
00025 srand( time(NULL) );
00026 int inSz = (rank+1)*( rand() % atoi(argv[1]) );
00027
00028 std::cout<<rank<<" inSz: "<<inSz<<std::endl;
00029
00030 Vec in;
00031 VecCreate(MPI_COMM_WORLD,&in);
00032 VecSetSizes(in,inSz,PETSC_DECIDE);
00033 std::cout<<rank<<" Created Vec (in) "<<std::endl;
00034 if(npes > 1) {
00035 VecSetType(in,VECMPI);
00036 std::cout<<rank<<" Type: MPI "<<std::endl;
00037 }else {
00038 VecSetType(in,VECSEQ);
00039 std::cout<<rank<<" Type: SEQ "<<std::endl;
00040 }
00041
00042 int offRank;
00043
00044 if(!rank) {
00045 srand( time(NULL) );
00046 offRank = ( rand() % npes );
00047 }
00048
00049 par::Mpi_Bcast<int>(&offRank,1,0,MPI_COMM_WORLD);
00050
00051 std::cout<<rank<<" offRank: "<<offRank<<std::endl;
00052
00053 int *inSizes = new int[npes];
00054 par::Mpi_Allgather<int>(&inSz, inSizes, 1, MPI_COMM_WORLD);
00055
00056 int *inDisps = new int[npes];
00057 inDisps[0] = 0;
00058 for(int i = 1; i < npes; i++) {
00059 inDisps[i] = inDisps[i-1] + inSizes[i-1];
00060 }
00061
00062
00063 int outSz = inSizes[(rank + offRank)%npes];
00064
00065 std::cout<<rank<<" outSz: "<<outSz<<std::endl;
00066
00067 Vec out;
00068 VecCreate(MPI_COMM_WORLD, &out);
00069 VecSetSizes(out,outSz,PETSC_DECIDE);
00070 std::cout<<rank<<" Created Vec (out) "<<std::endl;
00071 if(npes > 1) {
00072 VecSetType(out,VECMPI);
00073 std::cout<<rank<<" Type: MPI "<<std::endl;
00074 }else {
00075 VecSetType(out,VECSEQ);
00076 std::cout<<rank<<" Type: SEQ "<<std::endl;
00077 }
00078
00079 int *outSizes = new int[npes];
00080 par::Mpi_Allgather<int>(&outSz, outSizes, 1, MPI_COMM_WORLD);
00081
00082 int *outDisps = new int[npes];
00083 outDisps[0] = 0;
00084 for(int i = 1; i < npes; i++) {
00085 outDisps[i] = outDisps[i-1] + outSizes[i-1];
00086 }
00087
00088 DendroIntL totalInSz;
00089 DendroIntL totalOutSz;
00090
00091 DendroIntL inSzCopy = inSz;
00092 DendroIntL outSzCopy = outSz;
00093
00094 par::Mpi_Reduce<DendroIntL>(&inSzCopy, &totalInSz, 1, MPI_SUM,0, MPI_COMM_WORLD);
00095 par::Mpi_Reduce<DendroIntL>(&outSzCopy, &totalOutSz, 1, MPI_SUM,0, MPI_COMM_WORLD);
00096
00097 if(!rank) {
00098 assert(totalInSz == totalOutSz);
00099 std::cout<<"Global sizes of in and out match."<<std::endl;
00100 }
00101
00102
00103 PetscScalar * inArr;
00104 VecZeroEntries(in);
00105 VecGetArray(in,&inArr);
00106
00107 srand( time(NULL) );
00108 for(int i = 0; i < inSz; i++) {
00109 inArr[i] = ((double)(rank+1))*(((double)rand())/((double)RAND_MAX));
00110 }
00111
00112 PetscScalar *fullInArr = NULL;
00113 if(!rank) {
00114 fullInArr = new PetscScalar[totalInSz];
00115 }
00116
00117 MPI_Gatherv(inArr,inSz,par::Mpi_datatype<PetscScalar>::value(),
00118 fullInArr,inSizes,inDisps,
00119 par::Mpi_datatype<PetscScalar>::value(),0,MPI_COMM_WORLD);
00120
00121 VecRestoreArray(in,&inArr);
00122 delete [] inSizes;
00123 delete [] inDisps;
00124
00125 int *sendSz = NULL;
00126 int *sendOff = NULL;
00127 int *recvSz = NULL;
00128 int *recvOff = NULL;
00129
00130 ot::scatterValues(in, out, inSz, outSz, sendSz, sendOff, recvSz, recvOff, MPI_COMM_WORLD);
00131
00132 PetscInt tmpSzIn, tmpSzOut;
00133 VecGetLocalSize(in, &tmpSzIn);
00134 VecGetLocalSize(out, &tmpSzOut);
00135
00136 assert(tmpSzIn == inSz);
00137 assert(tmpSzOut == outSz);
00138
00139 MPI_Barrier(MPI_COMM_WORLD);
00140 if(!rank) {
00141 std::cout<<"The local sizes are preserved."<<std::endl;
00142 }
00143
00144 VecDestroy(in);
00145
00146 PetscScalar* outArr;
00147 VecGetArray(out,&outArr);
00148
00149 PetscScalar* fullOutArr = NULL;
00150 if(!rank) {
00151 fullOutArr = new PetscScalar[totalOutSz];
00152 }
00153
00154 MPI_Gatherv(outArr, outSz, par::Mpi_datatype<PetscScalar>::value(),
00155 fullOutArr, outSizes, outDisps, par::Mpi_datatype<PetscScalar>::value(),
00156 0,MPI_COMM_WORLD);
00157
00158 delete [] outSizes;
00159 delete [] outDisps;
00160 VecRestoreArray(out,&outArr);
00161 VecDestroy(out);
00162
00163 if(!rank) {
00164 for(DendroIntL i = 0; i < totalInSz; i++) {
00165 assert(fullInArr[i] == fullOutArr[i]);
00166 }
00167 delete [] fullInArr;
00168 delete [] fullOutArr;
00169 std::cout<<"Passed Vec Version."<<std::endl;
00170 }
00171 }
00172
00173 MPI_Barrier(MPI_COMM_WORLD);
00174
00175
00176 {
00177
00178 srand( time(NULL) );
00179 int inSz = (rank+1)*( rand() % atoi(argv[1]) );
00180
00181 std::cout<<rank<<" inSz: "<<inSz<<std::endl;
00182
00183 std::vector<int> in(inSz);
00184 std::cout<<rank<<" Created std::vec (in) "<<std::endl;
00185
00186 int offRank;
00187
00188 if(!rank) {
00189 srand( time(NULL) );
00190 offRank = ( rand() % npes );
00191 }
00192
00193 par::Mpi_Bcast<int>(&offRank,1,0,MPI_COMM_WORLD);
00194
00195 std::cout<<rank<<" offRank: "<<offRank<<std::endl;
00196
00197 int *inSizes = new int[npes];
00198 par::Mpi_Allgather<int>(&inSz, inSizes, 1, MPI_COMM_WORLD);
00199
00200 int *inDisps = new int[npes];
00201 inDisps[0] = 0;
00202 for(int i = 1; i < npes; i++) {
00203 inDisps[i] = inDisps[i-1] + inSizes[i-1];
00204 }
00205
00206
00207 int outSz = inSizes[(rank + offRank)%npes];
00208
00209 std::cout<<rank<<" outSz: "<<outSz<<std::endl;
00210
00211 std::vector<int> out(outSz);
00212 std::cout<<rank<<" Created std::vec (out) "<<std::endl;
00213
00214 int *outSizes = new int[npes];
00215 par::Mpi_Allgather<int>(&outSz, outSizes, 1, MPI_COMM_WORLD);
00216
00217 int *outDisps = new int[npes];
00218 outDisps[0] = 0;
00219 for(int i = 1; i < npes; i++) {
00220 outDisps[i] = outDisps[i-1] + outSizes[i-1];
00221 }
00222
00223 DendroIntL totalInSz;
00224 DendroIntL totalOutSz;
00225
00226 DendroIntL inSzCopy = inSz;
00227 DendroIntL outSzCopy = outSz;
00228
00229 par::Mpi_Reduce<DendroIntL>(&inSzCopy, &totalInSz, 1, MPI_SUM,0, MPI_COMM_WORLD);
00230 par::Mpi_Reduce<DendroIntL>(&outSzCopy, &totalOutSz, 1, MPI_SUM,0, MPI_COMM_WORLD);
00231
00232 if(!rank) {
00233 assert(totalInSz == totalOutSz);
00234 std::cout<<"Global sizes of in and out match."<<std::endl;
00235 }
00236
00237 srand( time(NULL) );
00238 for(int i = 0; i < inSz; i++) {
00239 in[i] = rand()*(rank+1);
00240 }
00241
00242 int *fullInArr = NULL;
00243 if(!rank) {
00244 fullInArr = new int[totalInSz];
00245 }
00246
00247 MPI_Gatherv(&(*in.begin()),inSz,MPI_INT,fullInArr,inSizes,inDisps,MPI_INT,0,MPI_COMM_WORLD);
00248
00249 delete [] inSizes;
00250 delete [] inDisps;
00251
00252 par::scatterValues<int>(in, out, outSz, MPI_COMM_WORLD);
00253
00254 int tmpSzIn = in.size();
00255 int tmpSzOut = out.size();
00256
00257 assert(tmpSzIn == inSz);
00258 assert(tmpSzOut == outSz);
00259
00260 MPI_Barrier(MPI_COMM_WORLD);
00261 if(!rank) {
00262 std::cout<<"The local sizes are preserved."<<std::endl;
00263 }
00264
00265 in.clear();
00266
00267 int *fullOutArr = NULL;
00268 if(!rank) {
00269 fullOutArr = new int[totalOutSz];
00270 }
00271
00272 MPI_Gatherv(&(*out.begin()),outSz,MPI_INT,fullOutArr,outSizes,outDisps,MPI_INT,0,MPI_COMM_WORLD);
00273
00274 delete [] outSizes;
00275 delete [] outDisps;
00276
00277 out.clear();
00278
00279 if(!rank) {
00280 for(int i = 0; i < totalInSz; i++) {
00281 assert(fullInArr[i] == fullOutArr[i]);
00282 }
00283 delete [] fullInArr;
00284 delete [] fullOutArr;
00285 std::cout<<"Passed STL Version."<<std::endl;
00286 }
00287 }
00288
00289 MPI_Barrier(MPI_COMM_WORLD);
00290
00291 PetscFinalize();
00292
00293 }
00294