Main Page | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members | Related Pages

compareMFiles.C

Go to the documentation of this file.
00001 
00002 #include "mpi.h"
00003 #include <iostream>
00004 #include "petsc.h"
00005 #include "sys.h"
00006 #include <vector>
00007 #include "parUtils.h"
00008 #include "testUtils.h"
00009 #include "externVars.h"
00010 #include "dendro.h"
00011 
00012 class MeshRecord {
00013         private:
00014                 unsigned int elementId;
00015                 unsigned int vtxId[8];
00016         public:
00017                 unsigned int getElementId() const {
00018                         return elementId;
00019                 }
00020 
00021                 unsigned int getVtxId(unsigned int vtxNum) const {
00022                         return vtxId[vtxNum];
00023                 }
00024 
00025                 MeshRecord() {
00026                         this->elementId = static_cast<unsigned int>(-1);
00027                         for(unsigned int j = 0; j < 8; j++) {
00028                                 this->vtxId[j] = static_cast<unsigned int>(-1);
00029                         }
00030                 }
00031 
00032                 MeshRecord(unsigned int eId, unsigned int* vId) {
00033                         this->elementId = eId;
00034                         for(unsigned int j = 0; j < 8; j++) {
00035                                 this->vtxId[j] = vId[j];
00036                         }
00037                 }
00038 
00039                 MeshRecord(MeshRecord const & other) {
00040                         this->elementId = other.elementId;
00041                         for(unsigned int j = 0; j < 8; j++) {
00042                                 ((this->vtxId)[j]) = ((other.vtxId)[j]);
00043                         }
00044                 }
00045 
00046                 MeshRecord & operator = (MeshRecord const  & other) {
00047                         if( this == (& other) ) {
00048                                 return (*this); 
00049                         }
00050                         this->elementId = other.elementId;
00051                         for(unsigned int j = 0; j < 8; j++) {
00052                                 ((this->vtxId)[j]) = ((other.vtxId)[j]);
00053                         }
00054                         return (*this);
00055                 }
00056 
00057                 bool  operator == ( MeshRecord const  &other) const {
00058                         bool result = (this->elementId == other.elementId);
00059                         for(unsigned int j = 0; ( (j < 8) && result ); j++) {
00060                                 result = (result && ( ((this->vtxId)[j]) == ((other.vtxId)[j]) ) );
00061                         }
00062                         return result;
00063                 }
00064 
00065                 bool  operator < ( MeshRecord const  &other) const {
00066                         //Sort using elementId first and then 0 thro 7 vtxId
00067                         //Since the elements should not repeat, eId should be unique.
00068                         if(this->elementId == other.elementId) {
00069                                 std::cout<<"WARNING: Comparing two records with the same elementId."<<std::endl;        
00070                                 unsigned int j = 0;
00071                                 while( ((this->vtxId)[j]) == ((other.vtxId)[j]) ) {
00072                                         j++;
00073                                 }
00074                                 return ( ((this->vtxId)[j]) < ((other.vtxId)[j]) );
00075                         } else {
00076                                 return (this->elementId < other.elementId);
00077                         }
00078                 }
00079 
00080                 bool  operator != (MeshRecord const  &other) const {
00081                         return (!((*this) == other));
00082                 }
00083 
00084                 bool  operator > ( MeshRecord const  &other) const {
00085                         return ( ((*this) != other) && ((*this) >= other) );
00086                 }
00087 
00088                 bool  operator <= ( MeshRecord const  &other) const {
00089                         return ( ((*this) < other) || ((*this) == other) );
00090                 }
00091 
00092                 bool  operator >= ( MeshRecord const  &other) const {
00093                         return (!((*this) < other));
00094                 }
00095 
00096                 friend std::ostream & operator << (std::ostream & os, MeshRecord const & v) ;
00097 }; //end class MeshRecord
00098 
00099 std::ostream & operator <<(std::ostream & os, MeshRecord const & r) {
00100         return (os << r.elementId <<" "<< r.vtxId[0] <<" "<<r.vtxId[1]<<" "
00101                         <<r.vtxId[2] <<" "<< r.vtxId[3] <<" "<<r.vtxId[4]<<" "
00102                         <<r.vtxId[5]<<" "<<r.vtxId[6]<<" "<<r.vtxId[7] );
00103 }//end fn.
00104 
00105 namespace par {
00106 
00107   template <typename T>
00108     class Mpi_datatype;
00109 
00110         template <>
00111                 class Mpi_datatype<MeshRecord> {
00112                         public: 
00113                                 static MPI_Datatype value() {
00114                                         static bool         first = true;
00115                                         static MPI_Datatype datatype;
00116                                         if (first) {
00117                                                 first = false;
00118                                                 MPI_Type_contiguous(sizeof(MeshRecord), MPI_BYTE, &datatype);
00119                                                 MPI_Type_commit(&datatype);
00120                                         }
00121                                         return datatype;
00122                                 }
00123                 };
00124 }//end namespace
00125 
00126 int main(int argc, char**argv) {
00127 
00128         PetscInitialize(&argc,&argv,0,0);
00129         ot::RegisterEvents();
00130 
00131         if(argc < 4) {
00132                 std::cout<<"exe fileBase file1Procs file2Procs "<<std::endl;
00133                 assert(false);
00134         }
00135 
00136         int rank, npes;
00137         MPI_Comm_rank(MPI_COMM_WORLD, &rank);
00138         MPI_Comm_size(MPI_COMM_WORLD, &npes);
00139 
00140         unsigned int file1Procs = atoi(argv[2]);
00141         unsigned int file2Procs = atoi(argv[3]);
00142 
00143         assert(file1Procs <= file2Procs);
00144         assert(file2Procs <= npes);
00145 
00146         std::vector<MeshRecord> vals1;
00147         std::vector<MeshRecord> vals2;
00148 
00149         if(rank < file1Procs) {
00150                 char fileName[100];
00151                 sprintf(fileName,"%s_%d_%d.out",argv[1],rank,file1Procs);
00152                 FILE * fptr = fopen(fileName,"rb");
00153                 assert(fptr != NULL);
00154                 unsigned int numElems = 0;
00155                 fread(&numElems, sizeof(unsigned int), 1, fptr);
00156                 vals1.resize(numElems);
00157                 for(unsigned int i = 0; i < numElems; i++ ) {   
00158                         unsigned int record[9];
00159                         fread(record, sizeof(unsigned int), 9, fptr);
00160                         vals1[i] = MeshRecord(record[0], (record+1));
00161                 }
00162                 fclose(fptr);
00163         }
00164 
00165         MPI_Barrier(MPI_COMM_WORLD);
00166 
00167         if(!rank) {
00168                 std::cout<<"Read File1 Successfully."<<std::endl;
00169         }
00170 
00171         MPI_Barrier(MPI_COMM_WORLD);
00172 
00173         if(rank < file2Procs) {
00174                 char fileName[100];
00175                 sprintf(fileName,"%s_%d_%d.out",argv[1],rank,file2Procs);
00176                 FILE * fptr = fopen(fileName,"rb");
00177                 assert(fptr != NULL);
00178                 unsigned int numElems = 0;
00179                 fread(&numElems, sizeof(unsigned int), 1, fptr);
00180                 vals2.resize(numElems);
00181                 for(unsigned int i = 0; i < numElems; i++ ) {   
00182                         unsigned int record[9];
00183                         fread(record, sizeof(unsigned int), 9, fptr);
00184                         vals2[i] = MeshRecord(record[0], (record+1));
00185                 }
00186                 fclose(fptr);
00187         }
00188 
00189         MPI_Barrier(MPI_COMM_WORLD);
00190 
00191         if(!rank) {
00192                 std::cout<<"Read File2 Successfully."<<std::endl;
00193         }
00194 
00195         MPI_Barrier(MPI_COMM_WORLD);
00196 
00197         DendroIntL locSz1 = vals1.size();
00198         DendroIntL globSz1 = 0;
00199         par::Mpi_Reduce<DendroIntL>(&locSz1, &globSz1, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00200 
00201         DendroIntL locSz2 = vals2.size();
00202         DendroIntL globSz2 = 0;
00203         par::Mpi_Reduce<DendroIntL>(&locSz2, &globSz2, 1, MPI_SUM, 0, MPI_COMM_WORLD);
00204 
00205         MPI_Barrier(MPI_COMM_WORLD);
00206 
00207         if(!rank) {
00208                 assert(globSz1 == globSz2);
00209                 std::cout << "Global Size: " << globSz1 << std::endl;
00210         }
00211 
00212         MPI_Barrier(MPI_COMM_WORLD);
00213 
00214         bool is1uniqueAndSorted =  par::test::isUniqueAndSorted<MeshRecord>(vals1, MPI_COMM_WORLD);
00215         bool is2uniqueAndSorted =  par::test::isUniqueAndSorted<MeshRecord>(vals2, MPI_COMM_WORLD);
00216 
00217         if(!rank) {
00218                 assert(is1uniqueAndSorted);
00219                 assert(is2uniqueAndSorted);
00220         }
00221 
00222         MPI_Barrier(MPI_COMM_WORLD);
00223 
00224         if(!rank) {
00225                 std::cout<<"The records are sorted and unique. "<<std::endl;
00226         }
00227 
00228         MPI_Barrier(MPI_COMM_WORLD);
00229 
00230         par::partitionW<MeshRecord>(vals1, NULL, MPI_COMM_WORLD);
00231         
00232         MPI_Barrier(MPI_COMM_WORLD);
00233 
00234         if(!rank) {
00235                 std::cout<<"Finished Partitioning List 1 "<<std::endl;
00236         }
00237 
00238         MPI_Barrier(MPI_COMM_WORLD);
00239 
00240         par::partitionW<MeshRecord>(vals2, NULL, MPI_COMM_WORLD);
00241 
00242         MPI_Barrier(MPI_COMM_WORLD);
00243 
00244         if(!rank) {
00245                 std::cout<<"Finished Partitioning List 2 "<<std::endl;
00246         }
00247 
00248         MPI_Barrier(MPI_COMM_WORLD);
00249 
00250         {
00251                 char outFile1[100];
00252                 sprintf(outFile1,"MeshList1_%d_%d.out",rank,npes);
00253                 FILE * fptr1 = fopen(outFile1,"wb");
00254                 unsigned int numElems = (vals1.size());
00255                 fwrite(&numElems, sizeof(unsigned int), 1, fptr1);
00256                 for(unsigned int i = 0; i < numElems; i++ ) {   
00257                         unsigned int record[9];
00258                         record[0] = vals1[i].getElementId();
00259                         for(unsigned int j = 0; j < 8; j++) {
00260                                 record[j+1] = vals1[i].getVtxId(j);
00261                         }
00262                         fwrite(record, sizeof(unsigned int), 9, fptr1);
00263                 }
00264                 fclose(fptr1);
00265         }
00266 
00267         {
00268                 char outFile2[100];
00269                 sprintf(outFile2,"MeshList2_%d_%d.out",rank,npes);
00270                 FILE * fptr2 = fopen(outFile2,"wb");
00271                 unsigned int numElems = (vals2.size());
00272                 fwrite(&numElems, sizeof(unsigned int), 1, fptr2);
00273                 for(unsigned int i = 0; i < numElems; i++ ) {   
00274                         unsigned int record[9];
00275                         record[0] = vals2[i].getElementId();
00276                         for(unsigned int j = 0; j < 8; j++) {
00277                                 record[j+1] = vals2[i].getVtxId(j);
00278                         }
00279                         fwrite(record, sizeof(unsigned int), 9, fptr2);
00280                 }
00281                 fclose(fptr2);
00282         }
00283 
00284         MPI_Barrier(MPI_COMM_WORLD);
00285 
00286         if(!rank) {
00287                 std::cout<<"Finished writing out the sorted records "<<std::endl;
00288         }
00289 
00290         MPI_Barrier(MPI_COMM_WORLD);
00291 
00292 
00293         assert(vals1.size() == vals2.size());
00294 
00295         MPI_Barrier(MPI_COMM_WORLD);
00296 
00297         bool Ipassed = false;
00298         if(vals1 == vals2) {
00299                 Ipassed = true;
00300         }
00301 
00302         bool *IpassedList = new bool[npes];
00303         par::Mpi_Allgather<bool>(&Ipassed, IpassedList, 1, MPI_COMM_WORLD);
00304 
00305         bool allPassed = true;
00306         unsigned int lastP = 0;
00307         for(unsigned int i = 0; ((i < npes) && allPassed); i++) {
00308                 allPassed = (allPassed && IpassedList[i]);
00309                 lastP = i;
00310         }
00311         delete [] IpassedList;
00312 
00313         if(!allPassed) {
00314                 if(rank == lastP) {
00315                         std::cout<<"First Failing on Proc: "<< lastP <<std::endl;
00316                         for(unsigned int i = 0; i < vals1.size(); i++) {
00317                                 if(vals1[i] != vals2[i]) {
00318                                         std::cout<<"i = "<<i<<" val1 = "<<vals1[i]<<" val2 = "<<vals2[i]<<std::endl;
00319                                         break;  
00320                                 }
00321                         }
00322                 }
00323         }
00324 
00325         vals1.clear();
00326         vals2.clear();
00327 
00328         MPI_Barrier(MPI_COMM_WORLD);
00329 
00330         if(!rank) {
00331                 if(allPassed) {
00332                         std::cout<<"The files match."<<std::endl;
00333                 }else {
00334                         std::cout<<"The files do not match."<<std::endl;
00335                 }
00336         }
00337 
00338         MPI_Barrier(MPI_COMM_WORLD);
00339 
00340         PetscFinalize();
00341 
00342 }
00343 

Generated on Tue Mar 24 16:14:07 2009 for DENDRO by  doxygen 1.3.9.1