00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014 #include <parmetislib.h>
00015
00016
00017 #define DEBUG_IPART_
00018
00019
00020
00021
00022
00023
00024
00025
00026 void InitMultisection(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00027 {
00028 int i, lpecut[2], gpecut[2], mypart, moptions[10];
00029 idxtype *vtxdist, *gwhere = NULL, *part, *label;
00030 GraphType *agraph;
00031 int *sendcounts, *displs;
00032 MPI_Comm newcomm, labelcomm;
00033
00034 IFSET(ctrl->dbglvl, DBG_TIME, starttimer(ctrl->InitPartTmr));
00035
00036
00037 agraph = AssembleMultisectedGraph(ctrl, graph, wspace);
00038 part = agraph->where;
00039 agraph->where = NULL;
00040
00041
00042 mypart = ctrl->mype%(ctrl->nparts/2);
00043 MPI_Comm_split(ctrl->comm, mypart, 0, &newcomm);
00044
00045
00046 agraph->ncon = 1;
00047 Moc_KeepPart(agraph, wspace, part, mypart);
00048 label = agraph->label;
00049 agraph->label = NULL;
00050
00051
00052 switch (ctrl->ipart) {
00053 case ISEP_EDGE:
00054 moptions[0] = 1;
00055 moptions[1] = 3;
00056 moptions[2] = 1;
00057 moptions[3] = 1;
00058 moptions[4] = 0;
00059 moptions[7] = ctrl->mype;
00060
00061 agraph->where = idxmalloc(agraph->nvtxs, "InitMultisection: agraph->where");
00062
00063 METIS_EdgeComputeSeparator(&agraph->nvtxs, agraph->xadj, agraph->adjncy,
00064 agraph->vwgt, agraph->adjwgt, moptions, &agraph->mincut, agraph->where);
00065 break;
00066 case ISEP_NODE:
00067 moptions[0] = 1;
00068 moptions[1] = 3;
00069 moptions[2] = 1;
00070 moptions[3] = 2;
00071 moptions[4] = 0;
00072 moptions[7] = ctrl->mype;
00073
00074 agraph->where = idxmalloc(agraph->nvtxs, "InitMultisection: agraph->where");
00075
00076 METIS_NodeComputeSeparator(&agraph->nvtxs, agraph->xadj, agraph->adjncy, agraph->vwgt,
00077 agraph->adjwgt, moptions, &agraph->mincut, agraph->where);
00078 break;
00079 default:
00080 errexit("Unknown ISEP type!\n");
00081 }
00082
00083 for (i=0; i<agraph->nvtxs; i++) {
00084 ASSERT(ctrl, agraph->where[i]>=0 && agraph->where[i]<=2);
00085 if (agraph->where[i] == 2)
00086 agraph->where[i] = ctrl->nparts+2*mypart;
00087 else
00088 agraph->where[i] += 2*mypart;
00089 }
00090
00091
00092 lpecut[0] = agraph->mincut;
00093 MPI_Comm_rank(newcomm, lpecut+1);
00094 MPI_Allreduce(lpecut, gpecut, 1, MPI_2INT, MPI_MINLOC, newcomm);
00095
00096
00097
00098
00099 if (lpecut[1] == gpecut[1] && gpecut[1] != 0)
00100 MPI_Send((void *)agraph->where, agraph->nvtxs, IDX_DATATYPE, 0, 1, newcomm);
00101 if (lpecut[1] == 0 && gpecut[1] != 0)
00102 MPI_Recv((void *)agraph->where, agraph->nvtxs, IDX_DATATYPE, gpecut[1], 1, newcomm, &ctrl->status);
00103
00104
00105 MPI_Comm_split(ctrl->comm, lpecut[1], 0, &labelcomm);
00106
00107
00108 if (lpecut[1] == 0) {
00109 gwhere = idxsmalloc(graph->gnvtxs, 0, "InitMultisection: gwhere");
00110 for (i=0; i<agraph->nvtxs; i++)
00111 gwhere[label[i]] = agraph->where[i];
00112 }
00113
00114 free(agraph->where);
00115 agraph->where = part;
00116
00117 if (lpecut[1] == 0) {
00118 MPI_Reduce((void *)gwhere, (void *)agraph->where, graph->gnvtxs, IDX_DATATYPE, MPI_SUM, 0, labelcomm);
00119 free(gwhere);
00120 }
00121
00122
00123 vtxdist = graph->vtxdist;
00124 ASSERT(ctrl, graph->where != NULL);
00125 free(graph->where);
00126 graph->where = idxmalloc(graph->nvtxs+graph->nrecv, "InitPartition: where");
00127
00128 sendcounts = imalloc(ctrl->npes, "InitPartitionNew: sendcounts");
00129 displs = imalloc(ctrl->npes, "InitPartitionNew: displs");
00130
00131 for (i=0; i<ctrl->npes; i++) {
00132 sendcounts[i] = vtxdist[i+1]-vtxdist[i];
00133 displs[i] = vtxdist[i];
00134 }
00135
00136 MPI_Scatterv((void *)agraph->where, sendcounts, displs, IDX_DATATYPE,
00137 (void *)graph->where, graph->nvtxs, IDX_DATATYPE, 0, ctrl->comm);
00138
00139 GKfree((void **)&sendcounts, (void **)&displs, (void **)&label, LTERM);
00140
00141 FreeGraph(agraph);
00142
00143 MPI_Comm_free(&newcomm);
00144 MPI_Comm_free(&labelcomm);
00145
00146 IFSET(ctrl->dbglvl, DBG_TIME, stoptimer(ctrl->InitPartTmr));
00147
00148 }
00149
00150
00151
00152
00153
00154
00155
00156 GraphType *AssembleMultisectedGraph(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00157 {
00158 int i, j, k, l, gnvtxs, nvtxs, gnedges, nedges, gsize;
00159 idxtype *xadj, *vwgt, *where, *adjncy, *adjwgt, *vtxdist, *imap;
00160 idxtype *axadj, *aadjncy, *aadjwgt, *avwgt, *awhere, *alabel;
00161 idxtype *mygraph, *ggraph;
00162 int *recvcounts, *displs, mysize;
00163 GraphType *agraph;
00164
00165 gnvtxs = graph->gnvtxs;
00166 nvtxs = graph->nvtxs;
00167 nedges = graph->xadj[nvtxs];
00168 xadj = graph->xadj;
00169 vwgt = graph->vwgt;
00170 where = graph->where;
00171 adjncy = graph->adjncy;
00172 adjwgt = graph->adjwgt;
00173 vtxdist = graph->vtxdist;
00174 imap = graph->imap;
00175
00176
00177 recvcounts = imalloc(ctrl->npes, "AssembleGraph: recvcounts");
00178 mysize = 3*nvtxs + 2*nedges;
00179 MPI_Allgather((void *)(&mysize), 1, MPI_INT, (void *)recvcounts, 1, MPI_INT, ctrl->comm);
00180
00181 displs = imalloc(ctrl->npes+1, "AssembleGraph: displs");
00182 displs[0] = 0;
00183 for (i=1; i<ctrl->npes+1; i++)
00184 displs[i] = displs[i-1] + recvcounts[i-1];
00185
00186
00187 mygraph = (mysize <= wspace->maxcore ? wspace->core : idxmalloc(mysize, "AssembleGraph: mygraph"));
00188 for (k=i=0; i<nvtxs; i++) {
00189 mygraph[k++] = xadj[i+1]-xadj[i];
00190 mygraph[k++] = vwgt[i];
00191 mygraph[k++] = where[i];
00192 for (j=xadj[i]; j<xadj[i+1]; j++) {
00193 mygraph[k++] = imap[adjncy[j]];
00194 mygraph[k++] = adjwgt[j];
00195 }
00196 }
00197 ASSERT(ctrl, mysize == k);
00198
00199
00200 gsize = displs[ctrl->npes];
00201 ggraph = (gsize <= wspace->maxcore-mysize ? wspace->core+mysize : idxmalloc(gsize, "AssembleGraph: ggraph"));
00202 MPI_Allgatherv((void *)mygraph, mysize, IDX_DATATYPE, (void *)ggraph, recvcounts, displs, IDX_DATATYPE, ctrl->comm);
00203
00204 GKfree((void **)&recvcounts, (void **)&displs, LTERM);
00205 if (mysize > wspace->maxcore)
00206 free(mygraph);
00207
00208 agraph = CreateGraph();
00209 agraph->nvtxs = gnvtxs;
00210 agraph->nedges = gnedges = (gsize-3*gnvtxs)/2;
00211
00212
00213 axadj = agraph->xadj = idxmalloc(gnvtxs+1, "AssembleGraph: axadj");
00214 avwgt = agraph->vwgt = idxmalloc(gnvtxs, "AssembleGraph: avwgt");
00215 awhere = agraph->where = idxmalloc(gnvtxs, "AssembleGraph: awhere");
00216 aadjncy = agraph->adjncy = idxmalloc(gnedges, "AssembleGraph: adjncy");
00217 aadjwgt = agraph->adjwgt = idxmalloc(gnedges, "AssembleGraph: adjwgt");
00218 alabel = agraph->label = idxmalloc(gnvtxs, "AssembleGraph: alabel");
00219
00220 for (k=j=i=0; i<gnvtxs; i++) {
00221 axadj[i] = ggraph[k++];
00222 avwgt[i] = ggraph[k++];
00223 awhere[i] = ggraph[k++];
00224 for (l=0; l<axadj[i]; l++) {
00225 aadjncy[j] = ggraph[k++];
00226 aadjwgt[j] = ggraph[k++];
00227 j++;
00228 }
00229 }
00230
00231
00232 MAKECSR(i, gnvtxs, axadj);
00233
00234 for (i=0; i<gnvtxs; i++)
00235 alabel[i] = i;
00236
00237 if (gsize > wspace->maxcore-mysize)
00238 free(ggraph);
00239
00240 return agraph;
00241 }
00242