Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
ElmerCSC
GitHub Repository: ElmerCSC/elmerfem
Path: blob/devel/elmergrid/src/metis-5.1.0/libmetis/pmetis.c
3206 views
1
/**
2
\file
3
\brief This file contains the top level routines for the multilevel recursive bisection
4
algorithm PMETIS.
5
6
\date Started 7/24/1997
7
\author George
8
\author Copyright 1997-2009, Regents of the University of Minnesota
9
\version\verbatim $Id: pmetis.c 10513 2011-07-07 22:06:03Z karypis $ \endverbatim
10
*/
11
12
13
#include "metislib.h"
14
15
16
/*************************************************************************/
17
/*! \ingroup api
18
\brief Recursive partitioning routine.
19
20
This function computes a partitioning of a graph based on multilevel
21
recursive bisection. It can be used to partition a graph into \e k
22
parts. The objective of the partitioning is to minimize the edgecut
23
subject to one or more balancing constraints.
24
25
\param[in] nvtxs is the number of vertices in the graph.
26
27
\param[in] ncon is the number of balancing constraints. For the standard
28
partitioning problem in which each vertex is either unweighted
29
or has a single weight, ncon should be 1.
30
31
\param[in] xadj is an array of size nvtxs+1 used to specify the starting
32
positions of the adjacency structure of the vertices in the
33
adjncy array.
34
35
\param[in] adjncy is an array of size to the sum of the degrees of the
36
graph that stores for each vertex the set of vertices that
37
is adjancent to.
38
39
\param[in] vwgt is an array of size nvtxs*ncon that stores the weights
40
of the vertices for each constraint. The ncon weights for the
41
ith vertex are stored in the ncon consecutive locations starting
42
at vwgt[i*ncon]. When ncon==1, a NULL value can be passed indicating
43
that all the vertices in the graph have the same weight.
44
45
\param[in] adjwgt is an array of size equal to adjncy, specifying the weight
46
for each edge (i.e., adjwgt[j] corresponds to the weight of the
47
edge stored in adjncy[j]).
48
A NULL value can be passed indicating that all the edges in the
49
graph have the same weight.
50
51
\param[in] nparts is the number of desired partitions.
52
53
\param[in] tpwgts is an array of size nparts*ncon that specifies the
54
desired weight for each part and constraint. The \e{target partition
55
weight} for the ith part and jth constraint is specified
56
at tpwgts[i*ncon+j] (the numbering of i and j starts from 0).
57
For each constraint, the sum of the tpwgts[] entries must be
58
1.0 (i.e., \f$ \sum_i tpwgts[i*ncon+j] = 1.0 \f$).
59
A NULL value can be passed indicating that the graph should
60
be equally divided among the parts.
61
62
\param[in] ubvec is an array of size ncon that specifies the allowed
63
load imbalance tolerance for each constraint.
64
For the ith part and jth constraint the allowed weight is the
65
ubvec[j]*tpwgts[i*ncon+j] fraction of the jth's constraint total
66
weight. The load imbalances must be greater than 1.0.
67
A NULL value can be passed indicating that the load imbalance
68
tolerance for each constraint should be 1.001 (for ncon==1)
69
or 1.01 (for ncon>1).
70
71
\params[in] options is the array for passing additional parameters
72
in order to customize the behaviour of the partitioning
73
algorithm.
74
75
\params[out] edgecut stores the cut of the partitioning.
76
77
\params[out] part is an array of size nvtxs used to store the
78
computed partitioning. The partition number for the ith
79
vertex is stored in part[i]. Based on the numflag parameter,
80
the numbering of the parts starts from either 0 or 1.
81
82
83
\returns
84
\retval METIS_OK indicates that the function returned normally.
85
\retval METIS_ERROR_INPUT indicates an input error.
86
\retval METIS_ERROR_MEMORY indicates that it could not allocate
87
the required memory.
88
89
*/
90
/*************************************************************************/
91
int METIS_PartGraphRecursive(idx_t *nvtxs, idx_t *ncon, idx_t *xadj,
92
idx_t *adjncy, idx_t *vwgt, idx_t *vsize, idx_t *adjwgt,
93
idx_t *nparts, real_t *tpwgts, real_t *ubvec, idx_t *options,
94
idx_t *objval, idx_t *part)
95
{
96
int sigrval=0, renumber=0;
97
graph_t *graph;
98
ctrl_t *ctrl;
99
100
/* set up malloc cleaning code and signal catchers */
101
if (!gk_malloc_init())
102
return METIS_ERROR_MEMORY;
103
104
gk_sigtrap();
105
106
if ((sigrval = gk_sigcatch()) != 0)
107
goto SIGTHROW;
108
109
110
/* set up the run parameters */
111
ctrl = SetupCtrl(METIS_OP_PMETIS, options, *ncon, *nparts, tpwgts, ubvec);
112
if (!ctrl) {
113
gk_siguntrap();
114
return METIS_ERROR_INPUT;
115
}
116
117
/* if required, change the numbering to 0 */
118
if (ctrl->numflag == 1) {
119
Change2CNumbering(*nvtxs, xadj, adjncy);
120
renumber = 1;
121
}
122
123
/* set up the graph */
124
graph = SetupGraph(ctrl, *nvtxs, *ncon, xadj, adjncy, vwgt, vsize, adjwgt);
125
126
/* allocate workspace memory */
127
AllocateWorkSpace(ctrl, graph);
128
129
/* start the partitioning */
130
IFSET(ctrl->dbglvl, METIS_DBG_TIME, InitTimers(ctrl));
131
IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_startcputimer(ctrl->TotalTmr));
132
133
*objval = MlevelRecursiveBisection(ctrl, graph, *nparts, part, ctrl->tpwgts, 0);
134
135
IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_stopcputimer(ctrl->TotalTmr));
136
IFSET(ctrl->dbglvl, METIS_DBG_TIME, PrintTimers(ctrl));
137
138
/* clean up */
139
FreeCtrl(&ctrl);
140
141
SIGTHROW:
142
/* if required, change the numbering back to 1 */
143
if (renumber)
144
Change2FNumbering(*nvtxs, xadj, adjncy, part);
145
146
gk_siguntrap();
147
gk_malloc_cleanup(0);
148
149
return metis_rcode(sigrval);
150
}
151
152
153
/*************************************************************************/
154
/*! This function is the top-level driver of the recursive bisection
155
routine. */
156
/*************************************************************************/
157
idx_t MlevelRecursiveBisection(ctrl_t *ctrl, graph_t *graph, idx_t nparts,
158
idx_t *part, real_t *tpwgts, idx_t fpart)
159
{
160
idx_t i, j, nvtxs, ncon, objval;
161
idx_t *label, *where;
162
graph_t *lgraph, *rgraph;
163
real_t wsum, *tpwgts2;
164
165
if ((nvtxs = graph->nvtxs) == 0) {
166
printf("\t***Cannot bisect a graph with 0 vertices!\n"
167
"\t***You are trying to partition a graph into too many parts!\n");
168
return 0;
169
}
170
171
ncon = graph->ncon;
172
173
/* determine the weights of the two partitions as a function of the weight of the
174
target partition weights */
175
WCOREPUSH;
176
tpwgts2 = rwspacemalloc(ctrl, 2*ncon);
177
for (i=0; i<ncon; i++) {
178
tpwgts2[i] = rsum((nparts>>1), tpwgts+i, ncon);
179
tpwgts2[ncon+i] = 1.0 - tpwgts2[i];
180
}
181
182
/* perform the bisection */
183
objval = MultilevelBisect(ctrl, graph, tpwgts2);
184
185
WCOREPOP;
186
187
label = graph->label;
188
where = graph->where;
189
for (i=0; i<nvtxs; i++)
190
part[label[i]] = where[i] + fpart;
191
192
if (nparts > 2)
193
SplitGraphPart(ctrl, graph, &lgraph, &rgraph);
194
195
/* Free the memory of the top level graph */
196
FreeGraph(&graph);
197
198
/* Scale the fractions in the tpwgts according to the true weight */
199
for (i=0; i<ncon; i++) {
200
wsum = rsum((nparts>>1), tpwgts+i, ncon);
201
rscale((nparts>>1), 1.0/wsum, tpwgts+i, ncon);
202
rscale(nparts-(nparts>>1), 1.0/(1.0-wsum), tpwgts+(nparts>>1)*ncon+i, ncon);
203
}
204
205
/* Do the recursive call */
206
if (nparts > 3) {
207
objval += MlevelRecursiveBisection(ctrl, lgraph, (nparts>>1), part,
208
tpwgts, fpart);
209
objval += MlevelRecursiveBisection(ctrl, rgraph, nparts-(nparts>>1), part,
210
tpwgts+(nparts>>1)*ncon, fpart+(nparts>>1));
211
}
212
else if (nparts == 3) {
213
FreeGraph(&lgraph);
214
objval += MlevelRecursiveBisection(ctrl, rgraph, nparts-(nparts>>1), part,
215
tpwgts+(nparts>>1)*ncon, fpart+(nparts>>1));
216
}
217
218
219
return objval;
220
}
221
222
223
/*************************************************************************/
224
/*! This function performs a multilevel bisection */
225
/*************************************************************************/
226
idx_t MultilevelBisect(ctrl_t *ctrl, graph_t *graph, real_t *tpwgts)
227
{
228
idx_t i, niparts, bestobj=0, curobj=0, *bestwhere=NULL;
229
graph_t *cgraph;
230
real_t bestbal=0.0, curbal=0.0;
231
232
Setup2WayBalMultipliers(ctrl, graph, tpwgts);
233
234
WCOREPUSH;
235
236
if (ctrl->ncuts > 1)
237
bestwhere = iwspacemalloc(ctrl, graph->nvtxs);
238
239
for (i=0; i<ctrl->ncuts; i++) {
240
cgraph = CoarsenGraph(ctrl, graph);
241
242
niparts = (cgraph->nvtxs <= ctrl->CoarsenTo ? SMALLNIPARTS : LARGENIPARTS);
243
Init2WayPartition(ctrl, cgraph, tpwgts, niparts);
244
245
Refine2Way(ctrl, graph, cgraph, tpwgts);
246
247
curobj = graph->mincut;
248
curbal = ComputeLoadImbalanceDiff(graph, 2, ctrl->pijbm, ctrl->ubfactors);
249
250
if (i == 0
251
|| (curbal <= 0.0005 && bestobj > curobj)
252
|| (bestbal > 0.0005 && curbal < bestbal)) {
253
bestobj = curobj;
254
bestbal = curbal;
255
if (i < ctrl->ncuts-1)
256
icopy(graph->nvtxs, graph->where, bestwhere);
257
}
258
259
if (bestobj == 0)
260
break;
261
262
if (i < ctrl->ncuts-1)
263
FreeRData(graph);
264
}
265
266
if (bestobj != curobj) {
267
icopy(graph->nvtxs, bestwhere, graph->where);
268
Compute2WayPartitionParams(ctrl, graph);
269
}
270
271
WCOREPOP;
272
273
return bestobj;
274
}
275
276
277
/*************************************************************************/
278
/*! This function splits a graph into two based on its bisection */
279
/*************************************************************************/
280
void SplitGraphPart(ctrl_t *ctrl, graph_t *graph, graph_t **r_lgraph,
281
graph_t **r_rgraph)
282
{
283
idx_t i, j, k, l, istart, iend, mypart, nvtxs, ncon, snvtxs[2], snedges[2];
284
idx_t *xadj, *vwgt, *adjncy, *adjwgt, *label, *where, *bndptr;
285
idx_t *sxadj[2], *svwgt[2], *sadjncy[2], *sadjwgt[2], *slabel[2];
286
idx_t *rename;
287
idx_t *auxadjncy, *auxadjwgt;
288
graph_t *lgraph, *rgraph;
289
290
WCOREPUSH;
291
292
IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_startcputimer(ctrl->SplitTmr));
293
294
nvtxs = graph->nvtxs;
295
ncon = graph->ncon;
296
xadj = graph->xadj;
297
vwgt = graph->vwgt;
298
adjncy = graph->adjncy;
299
adjwgt = graph->adjwgt;
300
label = graph->label;
301
where = graph->where;
302
bndptr = graph->bndptr;
303
304
ASSERT(bndptr != NULL);
305
306
rename = iwspacemalloc(ctrl, nvtxs);
307
308
snvtxs[0] = snvtxs[1] = snedges[0] = snedges[1] = 0;
309
for (i=0; i<nvtxs; i++) {
310
k = where[i];
311
rename[i] = snvtxs[k]++;
312
snedges[k] += xadj[i+1]-xadj[i];
313
}
314
315
lgraph = SetupSplitGraph(graph, snvtxs[0], snedges[0]);
316
sxadj[0] = lgraph->xadj;
317
svwgt[0] = lgraph->vwgt;
318
sadjncy[0] = lgraph->adjncy;
319
sadjwgt[0] = lgraph->adjwgt;
320
slabel[0] = lgraph->label;
321
322
rgraph = SetupSplitGraph(graph, snvtxs[1], snedges[1]);
323
sxadj[1] = rgraph->xadj;
324
svwgt[1] = rgraph->vwgt;
325
sadjncy[1] = rgraph->adjncy;
326
sadjwgt[1] = rgraph->adjwgt;
327
slabel[1] = rgraph->label;
328
329
snvtxs[0] = snvtxs[1] = snedges[0] = snedges[1] = 0;
330
sxadj[0][0] = sxadj[1][0] = 0;
331
for (i=0; i<nvtxs; i++) {
332
mypart = where[i];
333
334
istart = xadj[i];
335
iend = xadj[i+1];
336
if (bndptr[i] == -1) { /* This is an interior vertex */
337
auxadjncy = sadjncy[mypart] + snedges[mypart] - istart;
338
auxadjwgt = sadjwgt[mypart] + snedges[mypart] - istart;
339
for(j=istart; j<iend; j++) {
340
auxadjncy[j] = adjncy[j];
341
auxadjwgt[j] = adjwgt[j];
342
}
343
snedges[mypart] += iend-istart;
344
}
345
else {
346
auxadjncy = sadjncy[mypart];
347
auxadjwgt = sadjwgt[mypart];
348
l = snedges[mypart];
349
for (j=istart; j<iend; j++) {
350
k = adjncy[j];
351
if (where[k] == mypart) {
352
auxadjncy[l] = k;
353
auxadjwgt[l++] = adjwgt[j];
354
}
355
}
356
snedges[mypart] = l;
357
}
358
359
/* copy vertex weights */
360
for (k=0; k<ncon; k++)
361
svwgt[mypart][snvtxs[mypart]*ncon+k] = vwgt[i*ncon+k];
362
363
slabel[mypart][snvtxs[mypart]] = label[i];
364
sxadj[mypart][++snvtxs[mypart]] = snedges[mypart];
365
}
366
367
for (mypart=0; mypart<2; mypart++) {
368
iend = sxadj[mypart][snvtxs[mypart]];
369
auxadjncy = sadjncy[mypart];
370
for (i=0; i<iend; i++)
371
auxadjncy[i] = rename[auxadjncy[i]];
372
}
373
374
lgraph->nedges = snedges[0];
375
rgraph->nedges = snedges[1];
376
377
SetupGraph_tvwgt(lgraph);
378
SetupGraph_tvwgt(rgraph);
379
380
IFSET(ctrl->dbglvl, METIS_DBG_TIME, gk_stopcputimer(ctrl->SplitTmr));
381
382
*r_lgraph = lgraph;
383
*r_rgraph = rgraph;
384
385
WCOREPOP;
386
}
387
388
389