git-svn-id: svn://svn.icms.temple.edu/lammps-ro/trunk@7370 f3b2605a-c512-4ea7-a41b-209d697bcdaa

This commit is contained in:
sjplimp 2011-12-15 17:45:23 +00:00
parent 59165f7114
commit ee2ab8eac3
3 changed files with 25 additions and 17 deletions

View File

@ -84,7 +84,7 @@ void PairTersoffZBLOMP::read_file(char *file)
int params_per_line = 21;
char **words = new char*[params_per_line+1];
delete [] params;
memory->sfree(params);
params = NULL;
nparams = 0;

View File

@ -214,6 +214,8 @@ void Comm::set_proc_grid()
pmap->cart_map(1,procgrid,ncores,coregrid,myloc,procneigh,grid2proc);
else if (mapflag == XYZ)
pmap->xyz_map(xyz,procgrid,ncores,coregrid,myloc,procneigh,grid2proc);
// printf("AAA %d: %d %d: %d %d: %d %d\n",);
} else if (gridflag == NUMA) {
pmap->numa_map(0,coregrid,myloc,procneigh,grid2proc);
@ -1371,13 +1373,13 @@ void Comm::set_processors(int narg, char **arg)
if (iarg+6 > narg) error->all(FLERR,"Illegal processors command");
gridflag = TWOLEVEL;
ncores = atoi(arg[2]);
if (strcmp(arg[3],"*") == 0) user_coregrid[0] = 0;
else user_coregrid[0] = atoi(arg[3]);
if (strcmp(arg[4],"*") == 0) user_coregrid[1] = 0;
else user_coregrid[1] = atoi(arg[4]);
if (strcmp(arg[5],"*") == 0) user_coregrid[2] = 0;
else user_coregrid[2] = atoi(arg[5]);
ncores = atoi(arg[iarg+2]);
if (strcmp(arg[iarg+3],"*") == 0) user_coregrid[0] = 0;
else user_coregrid[0] = atoi(arg[iarg+3]);
if (strcmp(arg[iarg+4],"*") == 0) user_coregrid[1] = 0;
else user_coregrid[1] = atoi(arg[iarg+4]);
if (strcmp(arg[iarg+5],"*") == 0) user_coregrid[2] = 0;
else user_coregrid[2] = atoi(arg[iarg+5]);
if (ncores <= 0 || user_coregrid[0] < 0 ||
user_coregrid[1] < 0 || user_coregrid[2] < 0)
@ -1392,7 +1394,7 @@ void Comm::set_processors(int narg, char **arg)
gridflag = CUSTOM;
delete [] customfile;
int n = strlen(arg[iarg+2]) + 1;
customfile = new char(n);
customfile = new char[n];
strcpy(customfile,arg[iarg+2]);
iarg += 1;
@ -1454,7 +1456,7 @@ void Comm::set_processors(int narg, char **arg)
if (iarg+2 > narg) error->all(FLERR,"Illegal processors command");
delete [] outfile;
int n = strlen(arg[iarg+1]) + 1;
outfile = new char(n);
outfile = new char[n];
strcpy(outfile,arg[iarg+1]);
iarg += 2;

View File

@ -243,6 +243,9 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
best_factors(numapossible,numafactors,numagrid,
nodegrid[0],nodegrid[1],nodegrid[2]);
memory->destroy(numafactors);
memory->destroy(nodefactors);
// assign a unique id to each node
node_id = 0;
@ -392,6 +395,7 @@ void ProcMap::cart_map(int reorder, int *procgrid, int ncores, int *coregrid,
/* ----------------------------------------------------------------------
map processors to 3d grid in XYZ order
called by onelevel
------------------------------------------------------------------------- */
void ProcMap::xyz_map(char *xyz, int *procgrid,
@ -441,6 +445,7 @@ void ProcMap::xyz_map(char *xyz, int *procgrid,
/* ----------------------------------------------------------------------
map processors to 3d grid in XYZ order
respect sub-grid of cores within each node
called by twolevel
------------------------------------------------------------------------- */
void ProcMap::xyz_map(char *xyz, int *procgrid, int ncores, int *coregrid,
@ -460,15 +465,15 @@ void ProcMap::xyz_map(char *xyz, int *procgrid, int ncores, int *coregrid,
inode = i/coregrid[0];
jnode = j/coregrid[1];
knode = k/coregrid[2];
icore = i - inode*icore;
jcore = j - jnode*jcore;
kcore = k - knode*kcore;
icore = i % coregrid[0];
jcore = j % coregrid[1];
kcore = k % coregrid[2];
if (xyz[0] == 'x' && xyz[1] == 'y' && xyz[2] == 'z')
if (xyz[0] == 'x' && xyz[1] == 'y' && xyz[2] == 'z') {
grid2proc[i][j][k] = ncores *
(knode*nodegrid[1]*nodegrid[0] + jnode*nodegrid[0] + inode) +
(kcore*coregrid[1]*coregrid[0] + jcore*coregrid[0] + icore);
else if (xyz[0] == 'x' && xyz[1] == 'z' && xyz[2] == 'y')
} else if (xyz[0] == 'x' && xyz[1] == 'z' && xyz[2] == 'y')
grid2proc[i][j][k] = ncores *
(jnode*nodegrid[2]*nodegrid[0] + knode*nodegrid[0] + inode) +
(jcore*coregrid[2]*coregrid[0] + kcore*coregrid[0] + icore);
@ -560,14 +565,14 @@ void ProcMap::numa_map(int reorder, int *numagrid,
myloc[1] = myloc[1] * numagrid[1] + y_offset;
myloc[2] = myloc[2] * numagrid[2] + z_offset;
// allgather of locations to fill grid2proc
// allgather of myloc into gridi to fill grid2proc
int nprocs;
MPI_Comm_size(world,&nprocs);
int **gridi;
memory->create(gridi,nprocs,3,"comm:gridi");
MPI_Allgather(&myloc,3,MPI_INT,gridi[0],3,MPI_INT,world);
MPI_Allgather(myloc,3,MPI_INT,gridi[0],3,MPI_INT,world);
for (int i = 0; i < nprocs; i++)
grid2proc[gridi[i][0]][gridi[i][1]][gridi[i][2]] = i;
memory->destroy(gridi);
@ -859,6 +864,7 @@ int ProcMap::best_factors(int npossible, int **factors, int *best,
area[1]/factors[m][0]/factors[m][2] +
area[2]/factors[m][1]/factors[m][2];
if (surf < bestsurf) {
bestsurf = surf;
best[0] = factors[m][0];
best[1] = factors[m][1];
best[2] = factors[m][2];