Skip to content

Commit

Permalink
Remove the unbalance after global sort by globalId
Browse files Browse the repository at this point in the history
  • Loading branch information
thilinarmtb committed Feb 11, 2019
1 parent 6454010 commit 929fa2f
Show file tree
Hide file tree
Showing 4 changed files with 31 additions and 3 deletions.
1 change: 0 additions & 1 deletion Makefile
Original file line number Diff line number Diff line change
@@ -1,5 +1,4 @@
MPI?=1
VALGRIND?=1
DEBUG?=0
UNDERSCORE?=1
CC?=mpicc
Expand Down
30 changes: 29 additions & 1 deletion src/genmap-rsb.c
Original file line number Diff line number Diff line change
Expand Up @@ -220,6 +220,33 @@ int GenmapFiedler(GenmapHandle h, GenmapComm c, int maxIter, int global) {
return iter;
}

void GenmapSplitByGlobalId(GenmapHandle h) {
GenmapLong start = GenmapGetLocalStartIndex(h);
GenmapLong nel = GenmapGetNGlobalElements(h);
GenmapLong id = GenmapCommRank(GenmapGetLocalComm(h));
GenmapLong np = GenmapCommSize(GenmapGetLocalComm(h));
GenmapInt lelt = GenmapGetNLocalElements(h);
GenmapElements elements = GenmapGetElements(h);

GenmapInt pNel = (GenmapInt)(nel / np);
GenmapInt nrem = (GenmapInt)(nel - pNel * np);
GenmapInt idCount = 0;
while(idCount * pNel + ((idCount < nrem) ? idCount : nrem) < start)
idCount++;

GenmapLong upLimit = idCount * pNel + ((idCount < nrem) ? idCount : nrem);
GenmapLong downLimit = start;
do {
GenmapInt end = upLimit - start < lelt ? (GenmapInt)(upLimit - start) : lelt;
GenmapInt i;
for(i = (GenmapInt)(downLimit - start); i < end; i++)
elements[i].proc = idCount - 1;
downLimit = upLimit;
idCount++;
upLimit = idCount * pNel + ((idCount < nrem) ? idCount : nrem);
} while(downLimit - start < lelt);
}

void GenmapSplitByMedian(GenmapHandle h, int *bin) {
GenmapLong start = GenmapGetLocalStartIndex(h);
GenmapLong nel = GenmapGetNGlobalElements(h);
Expand Down Expand Up @@ -331,7 +358,6 @@ void GenmapRSB(GenmapHandle h) {
} while(ipass < npass && iter == maxIter);

GenmapBinSort(h, 0, &buf0);

int bin;
GenmapSplitByMedian(h, &bin);
GenmapTransferToBins(h, 0, &buf0);
Expand All @@ -342,6 +368,8 @@ void GenmapRSB(GenmapHandle h) {

#if defined(GENMAP_PAUL)
GenmapBinSort(h, 1, &buf0);
GenmapSplitByGlobalId(h);
GenmapTransferToBins(h, 1, &buf0);
#endif
}

Expand Down
1 change: 1 addition & 0 deletions src/parRSB.c
Original file line number Diff line number Diff line change
Expand Up @@ -47,6 +47,7 @@ int parRSB_partMesh(long long *egl, long long *vl, int *negl,
}

GenmapSetNLocalElements(h, (GenmapInt)neglcon);
assert(neglcon > 0);
GenmapScan(h, GenmapGetGlobalComm(h));
GenmapSetNVertices(h, nve);

Expand Down
2 changes: 1 addition & 1 deletion tests/con/conReader.h
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,7 @@ int conRead(char *fname, struct con *c, MPI_Comm comm) {
MPI_File_read_all(fh, hdr, sizeof(hdr), MPI_BYTE, MPI_STATUS_IGNORE);
char ver[5];
int nelgt, nelgv, nv;
sscanf(hdr, "%s %d %d %d", &ver, &nelgt, &nelgv, &nv);
sscanf(hdr, "%s %d %d %d", &ver[0], &nelgt, &nelgv, &nv);

float byte_test;
MPI_File_read_all(fh, &byte_test, 1, MPI_FLOAT, MPI_STATUS_IGNORE);
Expand Down

0 comments on commit 929fa2f

Please sign in to comment.