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

parUtils.C

Go to the documentation of this file.
00001 
00008 #include "mpi.h"
00009 #include "binUtils.h"
00010 #include "dtypes.h"
00011 #include "parUtils.h"
00012 
00013 #ifdef __DEBUG__
00014 #ifndef __DEBUG_PAR__
00015 #define __DEBUG_PAR__
00016 #endif
00017 #endif
00018 
00019 namespace par {
00020 
00021   unsigned int splitCommBinary( MPI_Comm orig_comm, MPI_Comm *new_comm) {
00022     int npes, rank;
00023 
00024     MPI_Group  orig_group, new_group;
00025 
00026     MPI_Comm_size(orig_comm, &npes);
00027     MPI_Comm_rank(orig_comm, &rank);
00028 
00029     unsigned int splitterRank = binOp::getPrevHighestPowerOfTwo(npes);
00030 
00031     int *ranksAsc, *ranksDesc;
00032     //Determine sizes for the 2 groups 
00033     ranksAsc = new int[splitterRank];
00034     ranksDesc = new int[( npes - splitterRank)];
00035 
00036     int numAsc = 0;
00037     int numDesc = ( npes - splitterRank - 1);
00038 
00039     //This is the main mapping between old ranks and new ranks.
00040     for(int i=0; i<npes; i++) {
00041       if( static_cast<unsigned int>(i) < splitterRank) {
00042         ranksAsc[numAsc] = i;
00043         numAsc++;
00044       }else {
00045         ranksDesc[numDesc] = i;
00046         numDesc--;
00047       }
00048     }//end for i
00049 
00050     MPI_Comm_group(orig_comm, &orig_group);
00051 
00052     /* Divide tasks into two distinct groups based upon rank */
00053     if (static_cast<unsigned int>(rank) < splitterRank) {
00054       MPI_Group_incl(orig_group, splitterRank, ranksAsc, &new_group);
00055     }else {
00056       MPI_Group_incl(orig_group, (npes-splitterRank), ranksDesc, &new_group);
00057     }
00058 
00059     MPI_Comm_create(orig_comm, new_group, new_comm);
00060 
00061     delete [] ranksAsc;
00062     delete [] ranksDesc;
00063 
00064     return splitterRank;
00065   }//end function
00066 
00067   unsigned int splitCommBinaryNoFlip( MPI_Comm orig_comm, MPI_Comm *new_comm) {
00068     int npes, rank;
00069 
00070     MPI_Group  orig_group, new_group;
00071 
00072     MPI_Comm_size(orig_comm, &npes);
00073     MPI_Comm_rank(orig_comm, &rank);
00074 
00075     unsigned int splitterRank =  binOp::getPrevHighestPowerOfTwo(npes);
00076 
00077     int *ranksAsc, *ranksDesc;
00078     //Determine sizes for the 2 groups 
00079     ranksAsc = new int[splitterRank];
00080     ranksDesc = new int[( npes - splitterRank)];
00081 
00082     int numAsc = 0;
00083     int numDesc = 0; //( npes - splitterRank - 1);
00084 
00085     //This is the main mapping between old ranks and new ranks.
00086     for(int i = 0; i < npes; i++) {
00087       if(static_cast<unsigned int>(i) < splitterRank) {
00088         ranksAsc[numAsc] = i;
00089         numAsc++;
00090       }else {
00091         ranksDesc[numDesc] = i;
00092         numDesc++;
00093       }
00094     }//end for i
00095 
00096     MPI_Comm_group(orig_comm, &orig_group);
00097 
00098     /* Divide tasks into two distinct groups based upon rank */
00099     if (static_cast<unsigned int>(rank) < splitterRank) {
00100       MPI_Group_incl(orig_group, splitterRank, ranksAsc, &new_group);
00101     }else {
00102       MPI_Group_incl(orig_group, (npes-splitterRank), ranksDesc, &new_group);
00103     }
00104 
00105     MPI_Comm_create(orig_comm, new_group, new_comm);
00106 
00107     delete [] ranksAsc;
00108     delete [] ranksDesc;
00109 
00110     return splitterRank;
00111   }//end function
00112 
00113   //create Comm groups and remove empty processors...
00114   int splitComm2way(bool iAmEmpty, MPI_Comm * new_comm, MPI_Comm comm) {
00115 #ifdef __PROFILE_WITH_BARRIER__
00116     MPI_Barrier(comm);
00117 #endif
00118     PROF_SPLIT_COMM_2WAY_BEGIN
00119 
00120       MPI_Group  orig_group, new_group;
00121     int size;
00122     MPI_Comm_size(comm, &size);
00123 
00124     bool* isEmptyList = new bool[size];
00125     par::Mpi_Allgather<bool>(&iAmEmpty, isEmptyList, 1, comm);
00126 
00127     int numActive=0, numIdle=0;
00128     for(int i = 0; i < size; i++) {
00129       if(isEmptyList[i]) {
00130         numIdle++;
00131       }else {
00132         numActive++;
00133       }
00134     }//end for i
00135 
00136     int* ranksActive = new int[numActive];
00137     int* ranksIdle = new int[numIdle];
00138 
00139     numActive=0;
00140     numIdle=0;
00141     for(int i = 0; i < size; i++) {
00142       if(isEmptyList[i]) {
00143         ranksIdle[numIdle] = i;
00144         numIdle++;
00145       }else {
00146         ranksActive[numActive] = i;
00147         numActive++;
00148       }
00149     }//end for i
00150 
00151     delete [] isEmptyList;      
00152 
00153     /* Extract the original group handle */
00154     MPI_Comm_group(comm, &orig_group);
00155 
00156     /* Divide tasks into two distinct groups based upon rank */
00157     if (!iAmEmpty) {
00158       MPI_Group_incl(orig_group, numActive, ranksActive, &new_group);
00159     }else {
00160       MPI_Group_incl(orig_group, numIdle, ranksIdle, &new_group);
00161     }
00162 
00163     /* Create new communicator */
00164     MPI_Comm_create(comm, new_group, new_comm);
00165 
00166     delete [] ranksActive;
00167     delete [] ranksIdle;
00168 
00169     PROF_SPLIT_COMM_2WAY_END
00170   }//end function
00171 
00172   int splitCommUsingSplittingRank(int splittingRank, MPI_Comm* new_comm,
00173       MPI_Comm comm) {
00174 #ifdef __PROFILE_WITH_BARRIER__
00175     MPI_Barrier(comm);
00176 #endif
00177     PROF_SPLIT_COMM_BEGIN
00178 
00179       MPI_Group  orig_group, new_group;
00180     int size;
00181     int rank;
00182     MPI_Comm_rank(comm, &rank);
00183     MPI_Comm_size(comm, &size);
00184 
00185     int* ranksActive = new int[splittingRank];
00186     int* ranksIdle = new int[size - splittingRank];
00187 
00188     for(int i = 0; i < splittingRank; i++) {
00189       ranksActive[i] = i;
00190     }
00191 
00192     for(int i = splittingRank; i < size; i++) {
00193       ranksIdle[i - splittingRank] = i;
00194     }
00195 
00196     /* Extract the original group handle */
00197     MPI_Comm_group(comm, &orig_group);
00198 
00199     /* Divide tasks into two distinct groups based upon rank */
00200     if (rank < splittingRank) {
00201       MPI_Group_incl(orig_group, splittingRank, ranksActive, &new_group);
00202     }else {
00203       MPI_Group_incl(orig_group, (size - splittingRank), ranksIdle, &new_group);
00204     }
00205 
00206     /* Create new communicator */
00207     MPI_Comm_create(comm, new_group, new_comm);
00208 
00209     delete [] ranksActive;
00210     delete [] ranksIdle;
00211 
00212     PROF_SPLIT_COMM_END
00213   }//end function
00214 
00215   //create Comm groups and remove empty processors...
00216   int splitComm2way(const bool* isEmptyList, MPI_Comm * new_comm, MPI_Comm comm) {
00217       
00218     MPI_Group  orig_group, new_group;
00219     int size, rank;
00220     MPI_Comm_size(comm, &size);
00221     MPI_Comm_rank(comm, &rank);
00222 
00223     int numActive=0, numIdle=0;
00224     for(int i = 0; i < size; i++) {
00225       if(isEmptyList[i]) {
00226         numIdle++;
00227       }else {
00228         numActive++;
00229       }
00230     }//end for i
00231 
00232     int* ranksActive = new int[numActive];
00233     int* ranksIdle = new int[numIdle];
00234 
00235     numActive=0;
00236     numIdle=0;
00237     for(int i = 0; i < size; i++) {
00238       if(isEmptyList[i]) {
00239         ranksIdle[numIdle] = i;
00240         numIdle++;
00241       }else {
00242         ranksActive[numActive] = i;
00243         numActive++;
00244       }
00245     }//end for i
00246 
00247     /* Extract the original group handle */
00248     MPI_Comm_group(comm, &orig_group);
00249 
00250     /* Divide tasks into two distinct groups based upon rank */
00251     if (!isEmptyList[rank]) {
00252       MPI_Group_incl(orig_group, numActive, ranksActive, &new_group);
00253     }else {
00254       MPI_Group_incl(orig_group, numIdle, ranksIdle, &new_group);
00255     }
00256 
00257     /* Create new communicator */
00258     MPI_Comm_create(comm, new_group, new_comm);
00259 
00260     delete [] ranksActive;
00261     delete [] ranksIdle;
00262 
00263     return 0;
00264   }//end function
00265 
00266 
00267 }// end namespace
00268 

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