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 GraphType *Moc_SetUpGraph(CtrlType *ctrl, int ncon, idxtype *vtxdist, idxtype *xadj,
00023 idxtype *vwgt, idxtype *adjncy, idxtype *adjwgt, int *wgtflag)
00024 {
00025 int i, j;
00026 GraphType *graph;
00027 int ltvwgts[MAXNCON];
00028
00029 graph = CreateGraph();
00030 graph->level = 0;
00031 graph->gnvtxs = vtxdist[ctrl->npes];
00032 graph->nvtxs = vtxdist[ctrl->mype+1]-vtxdist[ctrl->mype];
00033 graph->ncon = ncon;
00034 graph->nedges = xadj[graph->nvtxs];
00035 graph->xadj = xadj;
00036 graph->vwgt = vwgt;
00037 graph->adjncy = adjncy;
00038 graph->adjwgt = adjwgt;
00039 graph->vtxdist = vtxdist;
00040
00041
00042 if (((*wgtflag)&2) == 0)
00043 graph->vwgt = idxsmalloc(graph->nvtxs*ncon, 1, "Par_KMetis: vwgt");
00044
00045 if (((*wgtflag)&1) == 0)
00046 graph->adjwgt = idxsmalloc(graph->nedges, 1, "Par_KMetis: adjwgt");
00047
00048
00049 for (j=0; j<ncon; j++)
00050 ltvwgts[j] = 0;
00051
00052 for (i=0; i<graph->nvtxs; i++)
00053 for (j=0; j<ncon; j++)
00054 ltvwgts[j] += graph->vwgt[i*ncon+j];
00055
00056 for (j=0; j<ncon; j++)
00057 ctrl->tvwgts[j] = GlobalSESum(ctrl, ltvwgts[j]);
00058
00059
00060 for (i=0; i<ncon; i++) {
00061
00062 if (ctrl->tvwgts[i] == 0) {
00063 rprintf(ctrl, "ERROR: sum weight for constraint %d is zero\n", i);
00064 MPI_Finalize();
00065 exit(-1);
00066 }
00067 }
00068
00069
00070 graph->nvwgt = fmalloc(graph->nvtxs*ncon, "graph->nvwgt");
00071 for (i=0; i<graph->nvtxs; i++) {
00072 for (j=0; j<ncon; j++)
00073 graph->nvwgt[i*ncon+j] = (floattype)(graph->vwgt[i*ncon+j]) / (floattype)(ctrl->tvwgts[j]);
00074 }
00075
00076 srand(ctrl->seed);
00077
00078 return graph;
00079 }
00080
00081
00082
00083
00084
00085 void SetUpCtrl(CtrlType *ctrl, int nparts, int dbglvl, MPI_Comm comm)
00086 {
00087
00088 MPI_Comm_dup(comm, &(ctrl->gcomm));
00089 MPI_Comm_rank(ctrl->gcomm, &ctrl->mype);
00090 MPI_Comm_size(ctrl->gcomm, &ctrl->npes);
00091
00092 ctrl->dbglvl = dbglvl;
00093 ctrl->nparts = nparts;
00094 ctrl->comm = ctrl->gcomm;
00095 ctrl->xyztype = XYZ_SPFILL;
00096
00097 srand(ctrl->mype);
00098 }
00099
00100
00101
00102
00103
00104 void ChangeNumbering(idxtype *vtxdist, idxtype *xadj, idxtype *adjncy, idxtype *part, int npes, int mype, int from)
00105 {
00106 int i, nvtxs, nedges;
00107
00108 if (from == 1) {
00109 for (i=0; i<npes+1; i++)
00110 vtxdist[i]--;
00111
00112 nvtxs = vtxdist[mype+1]-vtxdist[mype];
00113 for (i=0; i<nvtxs+1; i++)
00114 xadj[i]--;
00115
00116 nedges = xadj[nvtxs];
00117 for (i=0; i<nedges; i++)
00118 adjncy[i]--;
00119 }
00120 else {
00121 nvtxs = vtxdist[mype+1]-vtxdist[mype];
00122 nedges = xadj[nvtxs];
00123
00124 for (i=0; i<npes+1; i++)
00125 vtxdist[i]++;
00126
00127 for (i=0; i<nvtxs+1; i++)
00128 xadj[i]++;
00129
00130 for (i=0; i<nedges; i++)
00131 adjncy[i]++;
00132
00133 for (i=0; i<nvtxs; i++)
00134 part[i]++;
00135
00136 }
00137 }
00138
00139
00140
00141
00142
00143 void ChangeNumberingMesh(idxtype *elmdist, idxtype *elements, idxtype *xadj,
00144 idxtype *adjncy, idxtype *part, int npes, int mype,
00145 int elmntlen, int from)
00146 {
00147 int i, nelms, nedges;
00148
00149 if (from == 1) {
00150 for (i=0; i<npes+1; i++)
00151 elmdist[i]--;
00152
00153 for (i=0; i<elmntlen; i++)
00154 elements[i]--;
00155 }
00156 else {
00157 nelms = elmdist[mype+1]-elmdist[mype];
00158 nedges = xadj[nelms];
00159
00160 for (i=0; i<npes+1; i++)
00161 elmdist[i]++;
00162
00163 for (i=0; i<elmntlen; i++)
00164 elements[i]++;
00165
00166 for (i=0; i<nelms+1; i++)
00167 xadj[i]++;
00168
00169 for (i=0; i<nedges; i++)
00170 adjncy[i]++;
00171
00172 if (part != NULL)
00173 for (i=0; i<nelms; i++)
00174 part[i]++;
00175 }
00176 }
00177
00178
00179
00180
00181
00182 void ChangeNumberingMesh2(idxtype *elmdist, idxtype *eptr, idxtype *eind,
00183 idxtype *xadj, idxtype *adjncy, idxtype *part,
00184 int npes, int mype, int from)
00185 {
00186 int i, nelms;
00187
00188 nelms = elmdist[mype+1]-elmdist[mype];
00189
00190 if (from == 1) {
00191 for (i=0; i<npes+1; i++)
00192 elmdist[i]--;
00193
00194 for (i=0; i<nelms+1; i++)
00195 eptr[i]--;
00196
00197 for (i=0; i<eptr[nelms]; i++)
00198 eind[i]--;
00199 }
00200 else {
00201 for (i=0; i<npes+1; i++)
00202 elmdist[i]++;
00203
00204 for (i=0; i<nelms+1; i++)
00205 eptr[i]++;
00206
00207 for (i=0; i<eptr[nelms]; i++)
00208 eind[i]++;
00209
00210 for (i=0; i<nelms+1; i++)
00211 xadj[i]++;
00212
00213 for (i=0; i<xadj[nelms]; i++)
00214 adjncy[i]++;
00215
00216 if (part != NULL)
00217 for (i=0; i<nelms; i++)
00218 part[i]++;
00219 }
00220 }
00221
00222
00223
00224
00225
00226
00227
00228 void GraphRandomPermute(GraphType *graph)
00229 {
00230 int i, j, k, tmp;
00231
00232 for (i=0; i<graph->nvtxs; i++) {
00233 for (j=graph->xadj[i]; j<graph->xadj[i+1]; j++) {
00234 k = graph->xadj[i] + RandomInRange(graph->xadj[i+1]-graph->xadj[i]);
00235 SWAP(graph->adjncy[j], graph->adjncy[k], tmp);
00236 SWAP(graph->adjwgt[j], graph->adjwgt[k], tmp);
00237 }
00238 }
00239 }
00240
00241
00242
00243
00244
00245
00246 void ComputeMoveStatistics(CtrlType *ctrl, GraphType *graph, int *nmoved, int *maxin, int *maxout)
00247 {
00248 int i, j, nvtxs;
00249 idxtype *vwgt, *where;
00250 idxtype *lpvtxs, *gpvtxs;
00251
00252 nvtxs = graph->nvtxs;
00253 vwgt = graph->vwgt;
00254 where = graph->where;
00255
00256 lpvtxs = idxsmalloc(ctrl->nparts, 0, "ComputeMoveStatistics: lpvtxs");
00257 gpvtxs = idxsmalloc(ctrl->nparts, 0, "ComputeMoveStatistics: gpvtxs");
00258
00259 for (j=i=0; i<nvtxs; i++) {
00260 lpvtxs[where[i]]++;
00261 if (where[i] != ctrl->mype)
00262 j++;
00263 }
00264
00265
00266
00267 MPI_Allreduce((void *)lpvtxs, (void *)gpvtxs, ctrl->nparts, IDX_DATATYPE, MPI_SUM, ctrl->comm);
00268
00269 *nmoved = GlobalSESum(ctrl, j);
00270 *maxout = GlobalSEMax(ctrl, j);
00271 *maxin = GlobalSEMax(ctrl, gpvtxs[ctrl->mype]-(nvtxs-j));
00272
00273 GKfree((void **)&lpvtxs, (void **)&gpvtxs, LTERM);
00274 }