00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015 #include <parmetislib.h>
00016
00017
00018
00019
00020
00021
00022 void PartitionSmallGraph(CtrlType *ctrl, GraphType *graph, WorkSpaceType *wspace)
00023 {
00024 int i, h, ncon, nparts, npes, mype;
00025 int moptions[10];
00026 int mynumflag, mywgtflag, me;
00027 idxtype *mypart;
00028 int lpecut[2], gpecut[2];
00029 GraphType *agraph;
00030 int *sendcounts, *displs;
00031 floattype *mytpwgts, *gnpwgts, *lnpwgts;
00032
00033 ncon = graph->ncon;
00034 nparts = ctrl->nparts;
00035
00036 MPI_Comm_size(ctrl->comm, &npes);
00037 MPI_Comm_rank(ctrl->comm, &mype);
00038
00039 SetUp(ctrl, graph, wspace);
00040 graph->where = idxmalloc(graph->nvtxs+graph->nrecv, "PartitionSmallGraph: where");
00041 agraph = Moc_AssembleAdaptiveGraph(ctrl, graph, wspace);
00042 mypart = idxmalloc(agraph->nvtxs, "mypart");
00043
00044 moptions[0] = 0;
00045 moptions[7] = ctrl->sync + mype;
00046 mynumflag = 0;
00047 mywgtflag = 3;
00048 if (ncon == 1) {
00049 METIS_WPartGraphKway2(&agraph->nvtxs, agraph->xadj, agraph->adjncy, agraph->vwgt,
00050 agraph->adjwgt, &mywgtflag, &mynumflag, &nparts, ctrl->tpwgts, moptions,
00051 &graph->mincut, mypart);
00052 }
00053 else {
00054 mytpwgts = fmalloc(nparts, "mytpwgts");
00055 for (i=0; i<nparts; i++)
00056 mytpwgts[i] = ctrl->tpwgts[i*ncon];
00057
00058 METIS_mCPartGraphRecursive2(&agraph->nvtxs, &ncon, agraph->xadj, agraph->adjncy,
00059 agraph->vwgt, agraph->adjwgt, &mywgtflag, &mynumflag, &nparts, mytpwgts,
00060 moptions, &graph->mincut, mypart);
00061
00062 free(mytpwgts);
00063 }
00064
00065 lpecut[0] = graph->mincut;
00066 lpecut[1] = mype;
00067 MPI_Allreduce(lpecut, gpecut, 1, MPI_2INT, MPI_MINLOC, ctrl->comm);
00068 graph->mincut = gpecut[0];
00069
00070 if (lpecut[1] == gpecut[1] && gpecut[1] != 0)
00071 MPI_Send((void *)mypart, agraph->nvtxs, IDX_DATATYPE, 0, 1, ctrl->comm);
00072 if (lpecut[1] == 0 && gpecut[1] != 0)
00073 MPI_Recv((void *)mypart, agraph->nvtxs, IDX_DATATYPE, gpecut[1], 1, ctrl->comm, &ctrl->status);
00074
00075 sendcounts = imalloc(npes, "sendcounts");
00076 displs = imalloc(npes, "displs");
00077
00078 for (i=0; i<npes; i++) {
00079 sendcounts[i] = graph->vtxdist[i+1]-graph->vtxdist[i];
00080 displs[i] = graph->vtxdist[i];
00081 }
00082
00083 MPI_Scatterv((void *)mypart, sendcounts, displs, IDX_DATATYPE,
00084 (void *)graph->where, graph->nvtxs, IDX_DATATYPE, 0, ctrl->comm);
00085
00086 lnpwgts = graph->lnpwgts = fmalloc(nparts*ncon, "lnpwgts");
00087 gnpwgts = graph->gnpwgts = fmalloc(nparts*ncon, "gnpwgts");
00088 sset(nparts*ncon, 0, lnpwgts);
00089 for (i=0; i<graph->nvtxs; i++) {
00090 me = graph->where[i];
00091 for (h=0; h<ncon; h++)
00092 lnpwgts[me*ncon+h] += graph->nvwgt[i*ncon+h];
00093 }
00094 MPI_Allreduce((void *)lnpwgts, (void *)gnpwgts, nparts*ncon, MPI_DOUBLE, MPI_SUM, ctrl->comm);
00095 GKfree((void**)&mypart, (void**)&sendcounts, (void**)&displs, LTERM);
00096 FreeGraph(agraph);
00097
00098 return;
00099 }
00100
00101
00102
00103
00104
00105
00106 void CheckInputs(int partType, int npes, int dbglvl, int *wgtflag, int *iwgtflag,
00107 int *numflag, int *inumflag, int *ncon, int *incon, int *nparts,
00108 int *inparts, floattype *tpwgts, floattype **itpwgts, floattype *ubvec,
00109 floattype *iubvec, floattype *ipc2redist, floattype *iipc2redist, int *options,
00110 int *ioptions, idxtype *part, MPI_Comm *comm)
00111 {
00112 int i, j;
00113 int doweabort, doiabort = 0;
00114 floattype tsum, *myitpwgts;
00115 int mgcnums[5] = {-1, 2, 3, 4, 2};
00116
00117
00118 if (part == NULL) {
00119 doiabort = 1;
00120 IFSET(dbglvl, DBG_INFO, printf("ERROR: part array is set to NULL.\n"));
00121 }
00122
00123
00124
00125
00126 if (wgtflag == NULL) {
00127 *iwgtflag = 0;
00128 IFSET(dbglvl, DBG_INFO, printf("WARNING: wgtflag is NULL. Using a value of 0.\n"));
00129 }
00130 else {
00131 *iwgtflag = *wgtflag;
00132 }
00133
00134
00135
00136
00137 if (numflag == NULL) {
00138 *inumflag = 0;
00139 IFSET(dbglvl, DBG_INFO, printf("WARNING: numflag is NULL. Using a value of 0.\n"));
00140 }
00141 else {
00142 if (*numflag != 0 && *numflag != 1) {
00143 IFSET(dbglvl, DBG_INFO, printf("WARNING: bad value for numflag %d. Using a value of 0.\n", *numflag));
00144 *inumflag = 0;
00145 }
00146 else {
00147 *inumflag = *numflag;
00148 }
00149 }
00150
00151
00152
00153
00154 if (ncon == NULL) {
00155 *incon = 1;
00156 IFSET(dbglvl, DBG_INFO, printf("WARNING: ncon is NULL. Using a value of 1.\n"));
00157 }
00158 else {
00159 if (*ncon < 1 || *ncon > MAXNCON) {
00160 IFSET(dbglvl, DBG_INFO, printf("WARNING: bad value for ncon %d. Using a value of 1.\n", *ncon));
00161 *incon = 1;
00162 }
00163 else {
00164 *incon = *ncon;
00165 }
00166 }
00167
00168
00169
00170
00171 if (nparts == NULL) {
00172 *inparts = npes;
00173 IFSET(dbglvl, DBG_INFO, printf("WARNING: nparts is NULL. Using a value of %d.\n", npes));
00174 }
00175 else {
00176 if (*nparts < 1 || *nparts > MAX_NPARTS) {
00177 IFSET(dbglvl, DBG_INFO, printf("WARNING: bad value for nparts %d. Using a value of %d.\n", *nparts, npes));
00178 *inparts = npes;
00179 }
00180 else {
00181 *inparts = *nparts;
00182 }
00183 }
00184
00185
00186
00187
00188 myitpwgts = *itpwgts = fmalloc((*inparts)*(*incon), "CheckInputs: itpwgts");
00189 if (tpwgts == NULL) {
00190 sset((*inparts)*(*incon), 1.0/(floattype)(*inparts), myitpwgts);
00191 IFSET(dbglvl, DBG_INFO, printf("WARNING: tpwgts is NULL. Setting all array elements to %.3f.\n", 1.0/(floattype)(*inparts)));
00192 }
00193 else {
00194 for (i=0; i<*incon; i++) {
00195 tsum = 0.0;
00196 for (j=0; j<*inparts; j++) {
00197 tsum += tpwgts[j*(*incon)+i];
00198 }
00199 if (fabs(1.0-tsum) < SMALLFLOAT)
00200 tsum = 1.0;
00201 for (j=0; j<*inparts; j++)
00202 myitpwgts[j*(*incon)+i] = tpwgts[j*(*incon)+i] / tsum;
00203 }
00204 }
00205
00206
00207
00208
00209 if (ubvec == NULL) {
00210 sset(*incon, 1.05, iubvec);
00211 IFSET(dbglvl, DBG_INFO, printf("WARNING: ubvec is NULL. Setting all array elements to 1.05.\n"));
00212 }
00213 else {
00214 for (i=0; i<*incon; i++) {
00215 if (ubvec[i] < 1.0 || ubvec[i] > (floattype)(*inparts)) {
00216 iubvec[i] = 1.05;
00217 IFSET(dbglvl, DBG_INFO, printf("WARNING: bad value for ubvec[%d]: %.3f. Setting value to 1.05.[%d]\n", i, ubvec[i], *inparts));
00218 }
00219 else {
00220 iubvec[i] = ubvec[i];
00221 }
00222 }
00223 }
00224
00225
00226
00227
00228 if (partType == ADAPTIVE_PARTITION) {
00229 if (ipc2redist != NULL) {
00230 if (*ipc2redist < SMALLFLOAT || *ipc2redist > 1000000.0) {
00231 IFSET(dbglvl, DBG_INFO, printf("WARNING: bad value for ipc2redist %.3f. Using a value of 1000.\n", *ipc2redist));
00232 *iipc2redist = 1000.0;
00233 }
00234 else {
00235 *iipc2redist = *ipc2redist;
00236 }
00237 }
00238 else {
00239 IFSET(dbglvl, DBG_INFO, printf("WARNING: ipc2redist is NULL. Using a value of 1000.\n"));
00240 *iipc2redist = 1000.0;
00241 }
00242 }
00243
00244
00245
00246
00247 if (options == NULL) {
00248 ioptions[0] = 0;
00249 IFSET(dbglvl, DBG_INFO, printf("WARNING: options is NULL. Using defaults\n"));
00250 }
00251 else {
00252 ioptions[0] = options[0];
00253 ioptions[1] = options[1];
00254 ioptions[2] = options[2];
00255 if (partType == ADAPTIVE_PARTITION || partType == REFINE_PARTITION)
00256 ioptions[3] = options[3];
00257 }
00258
00259
00260
00261
00262 if (comm == NULL) {
00263 IFSET(dbglvl, DBG_INFO, printf("ERROR: comm is NULL. Aborting\n"));
00264 abort();
00265 }
00266 else {
00267 MPI_Allreduce((void *)&doiabort, (void *)&doweabort, 1, MPI_INT, MPI_MAX, *comm);
00268 if (doweabort > 0)
00269 abort();
00270 }
00271
00272
00273 }
00274
00275