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

This commit is contained in:
sjplimp
2012-06-06 22:47:51 +00:00
parent f46eb9dedb
commit ef9e700545
1408 changed files with 58053 additions and 57983 deletions

View File

@ -5,7 +5,7 @@
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
@ -40,8 +40,8 @@ ProcMap::ProcMap(LAMMPS *lmp) : Pointers(lmp) {}
------------------------------------------------------------------------- */
void ProcMap::onelevel_grid(int nprocs, int *user_procgrid, int *procgrid,
int otherflag, int other_style,
int *other_procgrid, int *other_coregrid)
int otherflag, int other_style,
int *other_procgrid, int *other_coregrid)
{
int **factors;
@ -56,8 +56,8 @@ void ProcMap::onelevel_grid(int nprocs, int *user_procgrid, int *procgrid,
if (domain->dimension == 2) npossible = cull_2d(npossible,factors,3);
npossible = cull_user(npossible,factors,3,user_procgrid);
if (otherflag) npossible = cull_other(npossible,factors,3,
other_style,other_procgrid,
other_coregrid);
other_style,other_procgrid,
other_coregrid);
// user/other constraints make failure possible
@ -78,15 +78,15 @@ void ProcMap::onelevel_grid(int nprocs, int *user_procgrid, int *procgrid,
------------------------------------------------------------------------- */
void ProcMap::twolevel_grid(int nprocs, int *user_procgrid, int *procgrid,
int ncores, int *user_coregrid, int *coregrid,
int otherflag, int other_style,
int *other_procgrid, int *other_coregrid)
int ncores, int *user_coregrid, int *coregrid,
int otherflag, int other_style,
int *other_procgrid, int *other_coregrid)
{
int **nfactors,**cfactors,**factors;
if (nprocs % ncores)
if (nprocs % ncores)
error->all(FLERR,"Processors twogrid requires proc count "
"be a multiple of core count");
"be a multiple of core count");
// nfactors = list of all possible 3 factors of node count
// constrain by 2d
@ -117,8 +117,8 @@ void ProcMap::twolevel_grid(int nprocs, int *user_procgrid, int *procgrid,
npossible = cull_user(npossible,factors,4,user_procgrid);
if (otherflag) npossible = cull_other(npossible,factors,4,
other_style,other_procgrid,
other_coregrid);
other_style,other_procgrid,
other_coregrid);
// user/other constraints make failure possible
@ -147,7 +147,7 @@ void ProcMap::twolevel_grid(int nprocs, int *user_procgrid, int *procgrid,
------------------------------------------------------------------------- */
void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
int *numagrid)
int *numagrid)
{
// hardwire this for now
@ -163,7 +163,7 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
MPI_Allgather(node_name,MPI_MAX_PROCESSOR_NAME,MPI_CHAR,node_names,
MPI_MAX_PROCESSOR_NAME,MPI_CHAR,world);
std::string node_string = std::string(node_name);
// get number of procs per node
// NOTE: could do this without STL map
@ -177,17 +177,17 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
}
procs_per_node = name_map.begin()->second;
procs_per_numa = procs_per_node / numa_nodes;
delete [] node_names;
// error if any of these conditions met
if (nprocs % procs_per_numa || // total procs not a multiple of node
user_procgrid[0] > 1 || // user specified grid > 1 in any dim
user_procgrid[1] > 1 ||
user_procgrid[2] > 1)
error->all(FLERR,"Could not create numa grid of processors");
// user settings for the factorization per numa node
// currently not user settable
// if user specifies 1 for a proc grid dimension,
@ -199,7 +199,7 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
if (user_procgrid[0] == 1) user_numagrid[0] = 1;
if (user_procgrid[1] == 1) user_numagrid[1] = 1;
if (user_procgrid[2] == 1) user_numagrid[2] = 1;
// initial factorization within NUMA node
int **numafactors;
@ -207,7 +207,7 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
memory->create(numafactors,numapossible,3,"procmap:numafactors");
numapossible = factor(procs_per_numa,numafactors);
if (domain->dimension == 2)
if (domain->dimension == 2)
numapossible = cull_2d(numapossible,numafactors,3);
numapossible = cull_user(numapossible,numafactors,3,user_numagrid);
@ -232,7 +232,7 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
memory->create(nodefactors,nodepossible,3,"procmap:nodefactors");
nodepossible = factor(node_count,nodefactors);
if (domain->dimension == 2)
if (domain->dimension == 2)
nodepossible = cull_2d(nodepossible,nodefactors,3);
nodepossible = cull_user(nodepossible,nodefactors,3,user_nodegrid);
@ -240,14 +240,14 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
error->all(FLERR,"Could not create numa grid of processors");
best_factors(nodepossible,nodefactors,nodegrid,
numagrid[0],numagrid[1],numagrid[2]);
numagrid[0],numagrid[1],numagrid[2]);
// repeat NUMA node factorization using subdomain sizes
// refines the factorization if the user specified the node layout
// NOTE: this will not re-enforce user-procgrid constraint will it?
best_factors(numapossible,numafactors,numagrid,
nodegrid[0],nodegrid[1],nodegrid[2]);
nodegrid[0],nodegrid[1],nodegrid[2]);
memory->destroy(numafactors);
memory->destroy(nodefactors);
@ -273,7 +273,7 @@ void ProcMap::numa_grid(int nprocs, int *user_procgrid, int *procgrid,
------------------------------------------------------------------------- */
void ProcMap::custom_grid(char *cfile, int nprocs,
int *user_procgrid, int *procgrid)
int *user_procgrid, int *procgrid)
{
FILE *fp;
char line[MAXLINE];
@ -286,7 +286,7 @@ void ProcMap::custom_grid(char *cfile, int nprocs,
if (fp == NULL) error->one(FLERR,"Cannot open custom file");
// skip header = blank and comment lines
char *ptr;
if (!fgets(line,MAXLINE,fp))
error->one(FLERR,"Unexpected end of custom file");
@ -294,7 +294,7 @@ void ProcMap::custom_grid(char *cfile, int nprocs,
if (ptr = strchr(line,'#')) *ptr = '\0';
if (strspn(line," \t\n\r") != strlen(line)) break;
if (!fgets(line,MAXLINE,fp))
error->one(FLERR,"Unexpected end of custom file");
error->one(FLERR,"Unexpected end of custom file");
}
}
@ -320,9 +320,9 @@ void ProcMap::custom_grid(char *cfile, int nprocs,
if (me == 0) {
for (int i = 0; i < nprocs; i++) {
if (!fgets(line,MAXLINE,fp))
error->one(FLERR,"Unexpected end of custom file");
error->one(FLERR,"Unexpected end of custom file");
sscanf(line,"%d %d %d %d",
&cmap[i][0],&cmap[i][1],&cmap[i][2],&cmap[i][3]);
&cmap[i][0],&cmap[i][1],&cmap[i][2],&cmap[i][3]);
}
fclose(fp);
}
@ -349,12 +349,12 @@ void ProcMap::custom_grid(char *cfile, int nprocs,
------------------------------------------------------------------------- */
void ProcMap::cart_map(int reorder, int *procgrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
int periods[3];
periods[0] = periods[1] = periods[2] = 1;
MPI_Comm cartesian;
MPI_Cart_create(world,3,procgrid,periods,reorder,&cartesian);
MPI_Cart_get(cartesian,3,procgrid,periods,myloc);
MPI_Cart_shift(cartesian,0,1,&procneigh[0][0],&procneigh[0][1]);
@ -366,10 +366,10 @@ void ProcMap::cart_map(int reorder, int *procgrid,
for (i = 0; i < procgrid[0]; i++)
for (j = 0; j < procgrid[1]; j++)
for (k = 0; k < procgrid[2]; k++) {
coords[0] = i; coords[1] = j; coords[2] = k;
MPI_Cart_rank(cartesian,coords,&grid2proc[i][j][k]);
coords[0] = i; coords[1] = j; coords[2] = k;
MPI_Cart_rank(cartesian,coords,&grid2proc[i][j][k]);
}
MPI_Comm_free(&cartesian);
}
@ -380,7 +380,7 @@ void ProcMap::cart_map(int reorder, int *procgrid,
------------------------------------------------------------------------- */
void ProcMap::cart_map(int reorder, int *procgrid, int ncores, int *coregrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
// setup NUMA params that numa_grid() sets up
@ -405,7 +405,7 @@ void ProcMap::cart_map(int reorder, int *procgrid, int ncores, int *coregrid,
------------------------------------------------------------------------- */
void ProcMap::xyz_map(char *xyz, int *procgrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
int me;
MPI_Comm_rank(world,&me);
@ -414,22 +414,22 @@ void ProcMap::xyz_map(char *xyz, int *procgrid,
for (i = 0; i < procgrid[0]; i++)
for (j = 0; j < procgrid[1]; j++)
for (k = 0; k < procgrid[2]; k++) {
if (xyz[0] == 'x' && xyz[1] == 'y' && xyz[2] == 'z')
grid2proc[i][j][k] = k*procgrid[1]*procgrid[0] + j*procgrid[0] + i;
else if (xyz[0] == 'x' && xyz[1] == 'z' && xyz[2] == 'y')
grid2proc[i][j][k] = j*procgrid[2]*procgrid[0] + k*procgrid[0] + i;
else if (xyz[0] == 'y' && xyz[1] == 'x' && xyz[2] == 'z')
grid2proc[i][j][k] = k*procgrid[0]*procgrid[1] + i*procgrid[1] + j;
else if (xyz[0] == 'y' && xyz[1] == 'z' && xyz[2] == 'x')
grid2proc[i][j][k] = i*procgrid[2]*procgrid[1] + k*procgrid[1] + j;
else if (xyz[0] == 'z' && xyz[1] == 'x' && xyz[2] == 'y')
grid2proc[i][j][k] = j*procgrid[0]*procgrid[2] + i*procgrid[2] + k;
else if (xyz[0] == 'z' && xyz[1] == 'y' && xyz[2] == 'x')
grid2proc[i][j][k] = i*procgrid[1]*procgrid[2] + j*procgrid[2] + k;
if (xyz[0] == 'x' && xyz[1] == 'y' && xyz[2] == 'z')
grid2proc[i][j][k] = k*procgrid[1]*procgrid[0] + j*procgrid[0] + i;
else if (xyz[0] == 'x' && xyz[1] == 'z' && xyz[2] == 'y')
grid2proc[i][j][k] = j*procgrid[2]*procgrid[0] + k*procgrid[0] + i;
else if (xyz[0] == 'y' && xyz[1] == 'x' && xyz[2] == 'z')
grid2proc[i][j][k] = k*procgrid[0]*procgrid[1] + i*procgrid[1] + j;
else if (xyz[0] == 'y' && xyz[1] == 'z' && xyz[2] == 'x')
grid2proc[i][j][k] = i*procgrid[2]*procgrid[1] + k*procgrid[1] + j;
else if (xyz[0] == 'z' && xyz[1] == 'x' && xyz[2] == 'y')
grid2proc[i][j][k] = j*procgrid[0]*procgrid[2] + i*procgrid[2] + k;
else if (xyz[0] == 'z' && xyz[1] == 'y' && xyz[2] == 'x')
grid2proc[i][j][k] = i*procgrid[1]*procgrid[2] + j*procgrid[2] + k;
if (grid2proc[i][j][k] == me) {
myloc[0] = i; myloc[1] = j, myloc[2] = k;
}
if (grid2proc[i][j][k] == me) {
myloc[0] = i; myloc[1] = j, myloc[2] = k;
}
}
// proc IDs of neighbors
@ -455,7 +455,7 @@ void ProcMap::xyz_map(char *xyz, int *procgrid,
------------------------------------------------------------------------- */
void ProcMap::xyz_map(char *xyz, int *procgrid, int ncores, int *coregrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
int me;
MPI_Comm_rank(world,&me);
@ -468,41 +468,41 @@ void ProcMap::xyz_map(char *xyz, int *procgrid, int ncores, int *coregrid,
for (i = 0; i < procgrid[0]; i++)
for (j = 0; j < procgrid[1]; j++)
for (k = 0; k < procgrid[2]; k++) {
inode = i/coregrid[0];
jnode = j/coregrid[1];
knode = k/coregrid[2];
icore = i % coregrid[0];
jcore = j % coregrid[1];
kcore = k % coregrid[2];
inode = i/coregrid[0];
jnode = j/coregrid[1];
knode = k/coregrid[2];
icore = i % coregrid[0];
jcore = j % coregrid[1];
kcore = k % coregrid[2];
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')
grid2proc[i][j][k] = ncores *
(jnode*nodegrid[2]*nodegrid[0] + knode*nodegrid[0] + inode) +
(jcore*coregrid[2]*coregrid[0] + kcore*coregrid[0] + icore);
else if (xyz[0] == 'y' && xyz[1] == 'x' && xyz[2] == 'z')
grid2proc[i][j][k] = ncores *
(knode*nodegrid[0]*nodegrid[1] + inode*nodegrid[1] + jnode) +
(kcore*coregrid[0]*coregrid[1] + icore*coregrid[1] + jcore);
else if (xyz[0] == 'y' && xyz[1] == 'z' && xyz[2] == 'x')
grid2proc[i][j][k] = ncores *
(inode*nodegrid[2]*nodegrid[1] + knode*nodegrid[1] + jnode) +
(icore*coregrid[2]*coregrid[1] + kcore*coregrid[1] + jcore);
else if (xyz[0] == 'z' && xyz[1] == 'x' && xyz[2] == 'y')
grid2proc[i][j][k] = ncores *
(jnode*nodegrid[0]*nodegrid[2] + inode*nodegrid[2] + knode) +
(jcore*coregrid[0]*coregrid[2] + icore*coregrid[2] + kcore);
else if (xyz[0] == 'z' && xyz[1] == 'y' && xyz[2] == 'x')
grid2proc[i][j][k] = ncores *
(inode*nodegrid[1]*nodegrid[2] + jnode*nodegrid[2] + knode) +
(icore*coregrid[1]*coregrid[2] + jcore*coregrid[2] + kcore);
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')
grid2proc[i][j][k] = ncores *
(jnode*nodegrid[2]*nodegrid[0] + knode*nodegrid[0] + inode) +
(jcore*coregrid[2]*coregrid[0] + kcore*coregrid[0] + icore);
else if (xyz[0] == 'y' && xyz[1] == 'x' && xyz[2] == 'z')
grid2proc[i][j][k] = ncores *
(knode*nodegrid[0]*nodegrid[1] + inode*nodegrid[1] + jnode) +
(kcore*coregrid[0]*coregrid[1] + icore*coregrid[1] + jcore);
else if (xyz[0] == 'y' && xyz[1] == 'z' && xyz[2] == 'x')
grid2proc[i][j][k] = ncores *
(inode*nodegrid[2]*nodegrid[1] + knode*nodegrid[1] + jnode) +
(icore*coregrid[2]*coregrid[1] + kcore*coregrid[1] + jcore);
else if (xyz[0] == 'z' && xyz[1] == 'x' && xyz[2] == 'y')
grid2proc[i][j][k] = ncores *
(jnode*nodegrid[0]*nodegrid[2] + inode*nodegrid[2] + knode) +
(jcore*coregrid[0]*coregrid[2] + icore*coregrid[2] + kcore);
else if (xyz[0] == 'z' && xyz[1] == 'y' && xyz[2] == 'x')
grid2proc[i][j][k] = ncores *
(inode*nodegrid[1]*nodegrid[2] + jnode*nodegrid[2] + knode) +
(icore*coregrid[1]*coregrid[2] + jcore*coregrid[2] + kcore);
if (grid2proc[i][j][k] == me) {
myloc[0] = i; myloc[1] = j, myloc[2] = k;
}
if (grid2proc[i][j][k] == me) {
myloc[0] = i; myloc[1] = j, myloc[2] = k;
}
}
// proc IDs of neighbors
@ -526,28 +526,28 @@ void ProcMap::xyz_map(char *xyz, int *procgrid, int ncores, int *coregrid,
------------------------------------------------------------------------- */
void ProcMap::numa_map(int reorder, int *numagrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
// setup a per node communicator and find rank within
MPI_Comm node_comm;
MPI_Comm_split(world,node_id,0,&node_comm);
MPI_Comm_split(world,node_id,0,&node_comm);
int node_rank;
MPI_Comm_rank(node_comm,&node_rank);
// setup a per numa communicator and find rank within
MPI_Comm numa_comm;
int local_numa = node_rank / procs_per_numa;
MPI_Comm_split(node_comm,local_numa,0,&numa_comm);
MPI_Comm_split(node_comm,local_numa,0,&numa_comm);
int numa_rank;
MPI_Comm_rank(numa_comm,&numa_rank);
// setup a communicator with the rank 0 procs from each numa node
MPI_Comm numa_leaders;
MPI_Comm_split(world,numa_rank,0,&numa_leaders);
// use the MPI Cartesian routines to map the nodes to the grid
int periods[3];
@ -557,11 +557,11 @@ void ProcMap::numa_map(int reorder, int *numagrid,
MPI_Cart_create(numa_leaders,3,nodegrid,periods,reorder,&cartesian);
MPI_Cart_get(cartesian,3,nodegrid,periods,myloc);
}
// broadcast numa node location in grid to other procs in numa node
MPI_Bcast(myloc,3,MPI_INT,0,numa_comm);
// compute my location within the node grid
int z_offset = numa_rank / (numagrid[0] * numagrid[1]);
@ -570,7 +570,7 @@ void ProcMap::numa_map(int reorder, int *numagrid,
myloc[0] = myloc[0] * numagrid[0] + x_offset;
myloc[1] = myloc[1] * numagrid[1] + y_offset;
myloc[2] = myloc[2] * numagrid[2] + z_offset;
// allgather of myloc into gridi to fill grid2proc
int nprocs;
@ -582,7 +582,7 @@ void ProcMap::numa_map(int reorder, int *numagrid,
for (int i = 0; i < nprocs; i++)
grid2proc[gridi[i][0]][gridi[i][1]][gridi[i][2]] = i;
memory->destroy(gridi);
// proc IDs of neighbors
int minus,plus;
@ -611,7 +611,7 @@ void ProcMap::numa_map(int reorder, int *numagrid,
------------------------------------------------------------------------- */
void ProcMap::custom_map(int *procgrid,
int *myloc, int procneigh[3][2], int ***grid2proc)
int *myloc, int procneigh[3][2], int ***grid2proc)
{
int me,nprocs;
MPI_Comm_rank(world,&me);
@ -670,9 +670,9 @@ void ProcMap::output(char *file, int *procgrid, int ***grid2proc)
for (int i = 0; i < procgrid[0]; i++)
for (int j = 0; j < procgrid[1]; j++)
for (int k = 0; k < procgrid[2]; k++)
if (grid2proc[i][j][k] == me) {
ime = i; jme = j; kme = k;
}
if (grid2proc[i][j][k] == me) {
ime = i; jme = j; kme = k;
}
// polled comm of grid mapping info from each proc to proc 0
@ -695,16 +695,16 @@ void ProcMap::output(char *file, int *procgrid, int ***grid2proc)
if (me == 0) {
for (int iproc = 0; iproc < nprocs; iproc++) {
if (iproc) {
MPI_Send(&tmp,0,MPI_INT,iproc,0,world);
MPI_Recv(vec,6,MPI_INT,iproc,0,world,&status);
MPI_Recv(procname,MPI_MAX_PROCESSOR_NAME+1,MPI_CHAR,
iproc,0,world,&status);
MPI_Send(&tmp,0,MPI_INT,iproc,0,world);
MPI_Recv(vec,6,MPI_INT,iproc,0,world,&status);
MPI_Recv(procname,MPI_MAX_PROCESSOR_NAME+1,MPI_CHAR,
iproc,0,world,&status);
}
fprintf(fp,"%d %d %d: %d %d %d: %s\n",
vec[0],vec[1],vec[2],vec[3],vec[4],vec[5],procname);
vec[0],vec[1],vec[2],vec[3],vec[4],vec[5],procname);
}
} else {
MPI_Recv(&tmp,0,MPI_INT,0,0,world,&status);
MPI_Send(vec,6,MPI_INT,0,0,world);
@ -733,9 +733,9 @@ int ProcMap::factor(int n, int **factors)
for (j = 1; j <= nyz; j++) {
if (nyz % j) continue;
if (factors) {
factors[m][0] = i;
factors[m][1] = j;
factors[m][2] = nyz/j;
factors[m][0] = i;
factors[m][1] = j;
factors[m][2] = nyz/j;
}
m++;
}
@ -750,7 +750,7 @@ int ProcMap::factor(int n, int **factors)
------------------------------------------------------------------------- */
int ProcMap::combine_factors(int n1, int **factors1, int n2, int **factors2,
int **factors)
int **factors)
{
int m = 0;
for (int i = 0; i < n1; i++)
@ -806,9 +806,9 @@ int ProcMap::cull_user(int n, int **factors, int m, int *user_factors)
where Nx,Ny,Nz = node grid = procgrid/coregrid
------------------------------------------------------------------------- */
int ProcMap::cull_other(int n, int **factors, int m,
int other_style, int *other_procgrid,
int *other_coregrid)
int ProcMap::cull_other(int n, int **factors, int m,
int other_style, int *other_procgrid,
int *other_coregrid)
{
int i = 0;
while (i < n) {
@ -818,8 +818,8 @@ int ProcMap::cull_other(int n, int **factors, int m,
if ((other_procgrid[1]/other_coregrid[1]) % factors[i][1]) flag = 1;
if ((other_procgrid[2]/other_coregrid[2]) % factors[i][2]) flag = 1;
if (flag) {
for (int j = 0; j < m; j++) factors[i][j] = factors[n-1][j];
n--;
for (int j = 0; j < m; j++) factors[i][j] = factors[n-1][j];
n--;
} else i++;
}
}
@ -834,13 +834,13 @@ int ProcMap::cull_other(int n, int **factors, int m,
------------------------------------------------------------------------- */
int ProcMap::best_factors(int npossible, int **factors, int *best,
const int sx, const int sy, const int sz)
const int sx, const int sy, const int sz)
{
// determine cross-sectional areas for orthogonal and triclinic boxes
// for triclinic, area = cross product of 2 edge vectors stored in h matrix
// area[3] = surface area 3 box faces divided by sx,sy,sz
// area[0] = xy, area[1] = xz, area[2] = yz
double area[3];
if (domain->triclinic == 0) {
area[0] = domain->xprd * domain->yprd / (sx*sy);
@ -868,8 +868,8 @@ int ProcMap::best_factors(int npossible, int **factors, int *best,
double bestsurf = 2.0 * (area[0]+area[1]+area[2]);
for (int m = 0; m < npossible; m++) {
surf = area[0]/factors[m][0]/factors[m][1] +
area[1]/factors[m][0]/factors[m][2] +
surf = area[0]/factors[m][0]/factors[m][1] +
area[1]/factors[m][0]/factors[m][2] +
area[2]/factors[m][1]/factors[m][2];
if (surf < bestsurf) {
bestsurf = surf;
@ -887,11 +887,10 @@ int ProcMap::best_factors(int npossible, int **factors, int *best,
minus,plus = indices of neighboring processors in a dimension
------------------------------------------------------------------------- */
void ProcMap::grid_shift(int myloc, int nprocs, int &minus, int &plus)
void ProcMap::grid_shift(int myloc, int nprocs, int &minus, int &plus)
{
minus = myloc - 1;
if (minus < 0) minus = nprocs - 1;
plus = myloc + 1;
if (plus == nprocs) plus = 0;
}