Actual source code: lowstretch.cpp
1: #define PETSCKSP_DLL
3: #include <math.h>
4: #include <queue>
5: #include private/pcimpl.h
6: #include "boost/graph/adjacency_list.hpp"
7: #include "boost/graph/subgraph.hpp"
9: using namespace boost;
11: /*
12: Boost Graph type definitions
13: */
15: enum edge_length_t { edge_length };
16: enum edge_keep_t { edge_keep };
17: enum vertex_parent_t { vertex_parent };
18: enum vertex_children_t { vertex_children };
19: enum vertex_depth_t { vertex_depth };
20: namespace boost {
21: BOOST_INSTALL_PROPERTY(edge, length);
22: BOOST_INSTALL_PROPERTY(edge, keep);
23: BOOST_INSTALL_PROPERTY(vertex, parent);
24: BOOST_INSTALL_PROPERTY(vertex, children);
25: BOOST_INSTALL_PROPERTY(vertex, depth);
26: }
29: typedef property<vertex_parent_t, PetscInt,
30: property<vertex_children_t, std::vector<PetscInt>,
31: property<vertex_depth_t, PetscScalar,
32: property<vertex_index_t, PetscInt> > > > VertexProperty;
33: typedef property<edge_length_t, PetscScalar,
34: property<edge_keep_t, PetscTruth,
35: property<edge_index_t, PetscInt> > > EdgeProperty2;
36: typedef property<edge_weight_t, PetscScalar, EdgeProperty2> EdgeProperty;
37: // I wish I knew a better way to make a convenient edge property constructor
38: #define EDGE_PROPERTY(WEIGHT,LENGTH,KEEP) EdgeProperty(WEIGHT,EdgeProperty2(LENGTH,KEEP))
40: typedef subgraph<adjacency_list<vecS, vecS, undirectedS, VertexProperty, EdgeProperty> > Graph;
41: typedef graph_traits<Graph>::edge_descriptor Edge;
43: typedef property_map<Graph, edge_weight_t>::type EdgeWeight;
44: typedef property_map<Graph, edge_length_t>::type EdgeLength;
45: typedef property_map<Graph, edge_keep_t>::type EdgeKeep;
46: typedef property_map<Graph, edge_index_t>::type EdgeIndex;
47: typedef property_map<Graph, vertex_parent_t>::type VertexParent;
48: typedef property_map<Graph, vertex_children_t>::type VertexChildren;
49: typedef property_map<Graph, vertex_depth_t>::type VertexDepth;
51: typedef std::pair<PetscInt,PetscInt> PetscIntPair;
52: struct Component {
53: // static PetscInt next_id;
54: //static PetscInt max_id;
56: //PetscInt id;
57: std::vector<PetscInt> vertices; /* ordered s.t. root is first; parent precedes child */
58: /*
59: Component() {
60: id = next_id++;
61: }
62: */
64: };
65: struct ComponentPair {
66: Component *first;
67: Component *second;
68: std::vector<Edge> edges; // pointing from first to second
69: std::vector<std::pair<PetscScalar,PetscScalar> > lengthBelow;
70: std::pair<PetscScalar,PetscScalar> rootLengthBelow;
71: std::pair<PetscScalar,PetscScalar> rootCongestion;
73: ComponentPair() {
74: first = PETSC_NULL;
75: second = PETSC_NULL;
76: }
78: int getIndex(Component *c) {
79: if (first == c)
80: return 0;
81: else if (second == c)
82: return 1;
83: else
84: return -1;
85: }
87: Component* get(int i) {
88: return (i==0)?first:second;
89: }
91: void put(int i,Component *c) {
92: if (i==0)
93: first=c;
94: else
95: second=c;
96: }
98: bool match(Component *c1, Component *c2) {
99: return (first == c1 && second == c2) || (first == c2 && second == c1);
100: }
101: };
103: /* ShortestPathPriorityQueue is a priority queue in which each node (PQNode)
104: represents a potential shortest path to a vertex. Each node stores
105: the terminal vertex, the distance along the path, and (optionally)
106: the vertex (pred) adjacent to the terminal vertex.
107: The top node of the queue is the shortest path in the queue.
108: */
110: struct PQNode {
111: PetscInt vertex;
112: PetscInt pred;
113: PetscScalar dist;
115: PQNode() {}
117: PQNode(const PetscInt v,const PetscScalar d) {
118: vertex = v;
119: dist = d;
120: }
122: PQNode(const PetscInt v,const PetscInt p,const PetscScalar d) {
123: vertex = v;
124: pred = p;
125: dist = d;
126: }
128: bool operator<( const PQNode &a ) const {
129: return dist > a.dist;
130: }
131: };
132: typedef std::priority_queue<PQNode> ShortestPathPriorityQueue;
134: /*
135: Function headers
136: */
137: PetscErrorCode LowStretchSpanningTree(Mat mat,Mat *pre,
138: PetscScalar tol,PetscScalar& maxCong);
139: PetscErrorCode AugmentedLowStretchSpanningTree(Mat mat,Mat *pre,PetscTruth augment,
140: PetscScalar tol,PetscScalar& maxCong);
141: PetscErrorCode LowStretchSpanningTreeHelper(Graph& g,const PetscInt root,const PetscScalar alpha,PetscInt perm[]);
142: PetscErrorCode StarDecomp(const Graph g,const PetscInt root,const PetscScalar delta,const PetscScalar epsilon,
143: PetscInt& k,std::vector<PetscInt>& size,std::vector<std::vector<PetscInt> >& idx,
144: std::vector<PetscInt>& x,std::vector<PetscInt>& y);
145: PetscErrorCode AugmentSpanningTree(Graph& g,const PetscInt root,PetscScalar& maxCong);
146: PetscErrorCode DecomposeSpanningTree(Graph& g,const PetscInt root,
147: const PetscScalar maxInternalStretch,
148: const PetscScalar maxExternalWeight,
149: std::vector<Component*>& componentList,
150: std::vector<ComponentPair>& edgeComponentMap);
151: PetscErrorCode DecomposeSubTree(Graph& g,const PetscInt root,
152: const PetscScalar maxInternalStretch,
153: const PetscScalar maxExternalWeight,
154: std::vector<Component*>& componentList,
155: std::vector<ComponentPair>& edgeComponentMap,
156: Component*& currComponent,
157: PetscScalar& currInternalStretch,
158: PetscScalar& currExternalWeight);
159: PetscErrorCode AddBridges(Graph& g,
160: std::vector<Component*>& componentList,
161: std::vector<ComponentPair>& edgeComponentMap,
162: PetscScalar& maxCong);
165: /* -------------------------------------------------------------------------- */
166: /*
167: LowStretchSpanningTree - Applies EEST algorithm to construct
168: low-stretch spanning tree preconditioner
169:
171: Input Parameters:
172: . mat - input matrix
174: Output Parameter:
175: . pre - preconditioner matrix with cholesky factorization precomputed in place
176: */
179: PetscErrorCode LowStretchSpanningTree(Mat mat,Mat *prefact,
180: PetscScalar tol,PetscScalar& maxCong)
181: {
182: PetscErrorCode ierr;
186: AugmentedLowStretchSpanningTree(mat,prefact,PETSC_FALSE,tol,maxCong);
188: return(0);
189: }
190: /* -------------------------------------------------------------------------- */
191: /*
192: AugmentedLowStretchSpanningTree - Applies EEST algorithm to construct
193: low-stretch spanning tree preconditioner;
194: then augments tree with additional edges
195:
197: Input Parameters:
198: . mat - input matrix
199: . augment - augmenting options
201: Output Parameter:
202: . pre - preconditioner matrix with cholesky factorization precomputed in place
203: */
206: PetscErrorCode AugmentedLowStretchSpanningTree(Mat mat,Mat *prefact,PetscTruth augment,
207: PetscScalar tol,PetscScalar& maxCong)
208: {
209: PetscErrorCode ierr;
210: PetscInt *idx;
211: PetscInt start,end,n,ncols,i,j,k;
212: MatFactorInfo info;
213: // IS perm, iperm;
214: const PetscInt *cols_c;
215: const PetscScalar *vals_c;
216: PetscInt *rows, *cols, *dnz, *onz;
217: PetscScalar *vals, *diag, absval;
218: Mat *pre;
219: graph_traits<Graph>::out_edge_iterator e, e_end;
220: const PetscInt root = 0;
224: MatGetSize(mat,PETSC_NULL,&n);
226: Graph g(n);
228: EdgeKeep edge_keep_g = get(edge_keep_t(),g);
230: PetscMalloc3(n,PetscScalar,&diag,n,PetscInt,&dnz,n,PetscInt,&onz);
231: PetscMemzero(dnz, n * sizeof(PetscInt));
232: PetscMemzero(onz, n * sizeof(PetscInt));
233: MatGetOwnershipRange(mat, &start, &end);
234: for (i=0; i<n; i++) {
235: MatGetRow(mat,i,&ncols,&cols_c,&vals_c);
236: diag[i] = 0;
237: for (k=0; k<ncols; k++) {
238: if (cols_c[k] == i) {
239: diag[i] += vals_c[k];
240: } else if (vals_c[k] != 0) {
241: absval = vals_c[k]>0?vals_c[k]:-vals_c[k];
242: diag[i] -= absval;
243: if (cols_c[k] > i && absval > tol) {
244: // we set edge_weight = absval; edge_length = 1.0/absval; edge_keep = PETSC_FALSE
245: add_edge(i,cols_c[k],EDGE_PROPERTY(absval,1.0/absval,PETSC_FALSE),g);
246: }
247: }
248: }
249: MatRestoreRow(mat,i,&ncols,&cols_c,&vals_c);
250: }
252: PetscMalloc(n*sizeof(PetscInt),&idx);
253: put(get(vertex_depth_t(),g),root,0);
254: LowStretchSpanningTreeHelper(g,root,log(4.0/3)/(2.0*log(n)),idx);
256: if (augment) {
257: AugmentSpanningTree(g,root,maxCong);
258: }
260: pre = prefact;
261: MatCreate(PETSC_COMM_WORLD,pre);
262: MatSetSizes(*pre,PETSC_DECIDE,PETSC_DECIDE,n,n);
263: MatSetType(*pre,MATAIJ);
265: PetscMalloc3(1,PetscInt,&rows,n,PetscInt,&cols,n,PetscScalar,&vals);
266: for (i=0; i<n; i++) {
267: for (tie(e, e_end) = out_edges(i,g); e != e_end; e++) {
268: if (get(edge_keep_g,*e)) {
269: const PetscInt col = target(*e,g);
271: if (col >= start && col < end) {
272: dnz[i]++;
273: } else {
274: onz[i]++;
275: }
276: }
277: }
278: }
279: MatSeqAIJSetPreallocation(*pre, 0, dnz);
280: MatMPIAIJSetPreallocation(*pre, 0, dnz, 0, onz);
281: for (i=0; i<n; i++) {
282: rows[0] = i;
283: k = 0;
284: for (tie(e, e_end) = out_edges(i,g); e != e_end; e++) {
285: if (get(edge_keep_g,*e)) {
286: cols[k++] = target(*e,g);
287: }
288: }
289: MatGetValues(mat,1,rows,k,cols,vals);
290: for (j=0; j<k; j++) {
291: absval = vals[j]>0?vals[j]:-vals[j];
292: diag[i] += absval;
293: }
294: cols[k] = i;
295: vals[k] = diag[i];
296: MatSetValues(*pre,1,rows,k+1,cols,vals,INSERT_VALUES);
297: }
298: PetscFree3(rows,cols,vals);
299: PetscFree3(diag,dnz,onz);
301: MatAssemblyBegin(*pre,MAT_FINAL_ASSEMBLY);
302: MatAssemblyEnd(*pre,MAT_FINAL_ASSEMBLY);
304: /*
305: ISCreateGeneral(PETSC_COMM_WORLD,n,idx,&perm);
306: ISSetPermutation(perm);
307: ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
308: PetscFree(idx);
309: ISView(perm,PETSC_VIEWER_STDOUT_SELF);
310: */
311: IS rperm, cperm;
312: MatGetOrdering(*pre,MATORDERING_QMD,&rperm,&cperm);
314: MatFactorInfoInitialize(&info);
315: /*
316: MatCholeskyFactorSymbolic(*pre,iperm,&info,prefact);
317: MatCholeskyFactorNumeric(*pre,&info,prefact);
318: MatDestroy(*pre);
319: ISDestroy(perm);
320: ISDestroy(iperm);
321: */
322: MatLUFactor(*pre,rperm,cperm,&info);
323: ISDestroy(rperm);
324: ISDestroy(cperm);
326: return(0);
327: }
329: /* -------------------------------------------------------------------------- */
330: /*
331: LowStretchSpanningTreeHelper
333: Input Parameters:
334: . g - input graph; all edges have edge_keep = PETSC_FALSE;
335: root has vertex_depth set to distance from global root
336: . root - root vertex
337: . alpha - parameter
338: . perm - preallocated array of size num_vertices(g) in which to store vertex ordering
340: Output Parameter:
341: . g - edges in low-stretch spanning tree are marked with edge_keep = PETSC_TRUE;
342: also vertex_parent and vertex_children are set (vertex_parent of global root is undefined)
343: and vertex_depth is set to be distance from global root (where weight on edge is inverse distance)
344: . perm - list of vertices in which a vertex precedes its parent (with vertices referred to by global index)
345: */
348: PetscErrorCode LowStretchSpanningTreeHelper(Graph& g,const PetscInt root,const PetscScalar alpha,PetscInt perm[])
349: {
350: PetscInt n,i,j,k;
351: std::vector<PetscInt> size,x,y;
352: std::vector<std::vector<PetscInt> > idx;
357: EdgeLength edge_length_g = get(edge_length_t(),g);
358: EdgeKeep edge_keep_g = get(edge_keep_t(),g);
359: VertexParent vertex_parent_g = get(vertex_parent_t(),g);
360: VertexChildren vertex_children_g = get(vertex_children_t(),g);
361: VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
362: n = num_vertices(g);
364: if (n > 2) {
365: StarDecomp(g,root,1.0/3,alpha,k,size,idx,x,y);
366: j = n - size[0];
367: Graph& g1 = g.create_subgraph(idx[0].begin(),idx[0].end());
368: LowStretchSpanningTreeHelper(g1,g1.global_to_local(g.local_to_global(root)),alpha,perm+j);
369: for (i=1;i<=k;i++) {
370: Edge e = edge(x[i-1],y[i-1],g).first;
371: put(edge_keep_g,e,PETSC_TRUE);
372: put(vertex_parent_g,x[i-1],g.local_to_global(y[i-1]));
373: get(vertex_children_g,y[i-1]).push_back(g.local_to_global(x[i-1]));
374: put(vertex_depth_g,x[i-1],get(vertex_depth_g,y[i-1])+get(edge_length_g,e));
376: j -= size[i];
377: Graph& g1 = g.create_subgraph(idx[i].begin(),idx[i].end());
378: LowStretchSpanningTreeHelper(g1,g1.global_to_local(g.local_to_global(x[i-1])),alpha,perm+j);
379: }
380: } else if (n == 2) {
381: Edge e = *(out_edges(root,g).first);
382: PetscInt t = target(e,g);
384: put(edge_keep_g,e,PETSC_TRUE);
385: put(vertex_parent_g,t,g.local_to_global(root));
386: get(vertex_children_g,root).push_back(g.local_to_global(t));
387: put(vertex_depth_g,t,get(vertex_depth_g,root)+get(edge_length_g,e));
389: perm[0] = g.local_to_global(t);
390: perm[1] = g.local_to_global(root);
391: } else /* n == 1 */ {
392: perm[0] = g.local_to_global(root);
393: }
395: return(0);
396: }
400: /* -------------------------------------------------------------------------- */
401: /*
402: StarDecomp - calculate a star decomposition of the graph
404: Input Parameters:
405: . g - input graph
406: . root - center of star decomposition
407: . delta, epsilon - parameters of star decomposition
409: Output Parameter:
410: . k - number of partitions, not-including central one.
411: . size[i] - number of vertices in partition i
412: . idx[i] - list of vertices in partition i; partition 0 contains root
413: . x[i-1] - "anchor" vertex of non-central partition i
414: . y[i-i] - vertex in partition 0 forming "bridge" with x[i-1]
415: */
418: PetscErrorCode StarDecomp(Graph g,const PetscInt root,const PetscScalar delta,const PetscScalar epsilon,
419: PetscInt& k,std::vector<PetscInt>& size,std::vector<std::vector<PetscInt> >& idx,
420: std::vector<PetscInt>& x,std::vector<PetscInt>& y)
421: {
422: PetscInt n,m,edgesLeft;
424: ShortestPathPriorityQueue pq;
425: PetscScalar radius;
426: PetscInt centerSize;
427: std::vector<PetscInt> centerIdx;
428: PQNode node;
432: EdgeWeight edge_weight_g = get(edge_weight_t(),g);
433: EdgeLength edge_length_g = get(edge_length_t(),g);
434: n = num_vertices(g);
435: m = num_edges(g);
436: edgesLeft = m;
438: std::vector<PetscInt> pred(n,-1);
439: std::vector<PetscInt> succ[n];
440: std::vector<PetscInt>::iterator i;
441: PetscScalar dist[n];
442: std::vector<PetscTruth> taken(n,PETSC_FALSE);
444: /** form tree of shortest paths to root **/
445: graph_traits<Graph>::out_edge_iterator e, e_end;
446: for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
447: PetscInt t = target(*e,g);
448: pq.push(PQNode(t,root,get(edge_length_g,*e)));
449: }
450: pred[root] = root;
451: while (!pq.empty()) {
452: node = pq.top();pq.pop();
453: if (pred[node.vertex] == -1) {
454: succ[node.pred].push_back(node.vertex);
455: pred[node.vertex] = node.pred;
456: dist[node.vertex] = node.dist;
457: for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
458: PetscInt t = target(*e,g);
459: if (pred[t] == -1) {
460: pq.push(PQNode(t,node.vertex,node.dist+get(edge_length_g,*e)));
461: }
462: }
463: radius = node.dist;
464: }
465: }
467: /** BALL CUT **/
468: for (i=succ[root].begin();i!=succ[root].end();i++) {
469: pq.push(PQNode(*i,dist[*i]));
470: }
471: PetscScalar boundary = 0;
472: PetscInt edgeCount = 0;
473: centerIdx.push_back(g.local_to_global(root));
474: taken[root] = PETSC_TRUE;
475: centerSize = 1;
476: for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
477: boundary += get(edge_weight_g,*e);
478: edgeCount++;
479: }
480: const PetscScalar minRadius = delta*radius;
481: while (dist[pq.top().vertex] < minRadius) {
482: assert(!pq.empty());
483: node = pq.top();pq.pop();
484: centerIdx.push_back(g.local_to_global(node.vertex));
485: taken[node.vertex] = PETSC_TRUE;
486: centerSize++;
487: for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
488: if (taken[target(*e,g)]) {
489: boundary -= get(edge_weight_g,*e);
490: } else {
491: boundary += get(edge_weight_g,*e);
492: edgeCount++;
493: }
494: }
495: for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
496: pq.push(PQNode(*i,dist[*i]));
497: }
498: }
499: while (boundary > (edgeCount+1)*log(m)/(log(2)*(1-2*delta)*radius)) {
500: assert(!pq.empty());
501: node = pq.top();pq.pop();
502: centerIdx.push_back(g.local_to_global(node.vertex));
503: taken[node.vertex] = PETSC_TRUE;
504: centerSize++;
505: for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
506: if (taken[target(*e,g)]) {
507: boundary -= get(edge_weight_g,*e);
508: } else {
509: boundary += get(edge_weight_g,*e);
510: edgeCount++;
511: }
512: }
513: for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
514: pq.push(PQNode(*i,dist[*i]));
515: }
516: }
517: size.push_back(centerSize);
518: idx.push_back(centerIdx);
519: edgesLeft -= edgeCount;
521: k = 0;
522: assert(!pq.empty());
523: std::queue<PetscInt> anchor_q;
524: ShortestPathPriorityQueue cone_pq;
525: std::vector<PetscInt> cone_succ[n];
526: std::vector<PetscTruth> cone_found(n,PETSC_FALSE);
528: /** form tree of shortest paths to an anchor **/
529: while (!pq.empty()) {
530: node = pq.top();pq.pop();
531: cone_found[node.vertex] = PETSC_TRUE;
532: anchor_q.push(node.vertex);
533: for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
534: PetscInt t = target(*e,g);
535: if (!taken[t]) {
536: cone_pq.push(PQNode(t,node.vertex,get(edge_length_g,*e)));
537: }
538: }
539: }
540: while (!cone_pq.empty()) {
541: node = cone_pq.top();cone_pq.pop();
542: if (!cone_found[node.vertex]) {
543: cone_succ[node.pred].push_back(node.vertex);
544: cone_found[node.vertex] = PETSC_TRUE;
545: for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
546: PetscInt t = target(*e,g);
547: if (!taken[t] && !cone_found[t]) {
548: cone_pq.push(PQNode(t,node.vertex,node.dist+get(edge_length_g,*e)));
549: }
550: }
551: }
552: }
554: while (!anchor_q.empty()) {
555: /** CONE CUT **/
556: PetscInt anchor = anchor_q.front();anchor_q.pop();
557: if (!taken[anchor]) {
558: PetscInt v;
559: PetscInt thisSize = 0;
560: std::vector<PetscInt> thisIdx;
561: std::queue<PetscInt> q;
562: ShortestPathPriorityQueue mycone_pq;
563: std::vector<PetscTruth> mycone_taken(n,PETSC_FALSE);
564: PetscInt initialInternalConeEdges = 0;
566: boundary = 0;
567: edgeCount = 0;
568: q.push(anchor);
569: while (!q.empty()) {
570: v = q.front();q.pop();
571: taken[v] = PETSC_TRUE;
572: mycone_taken[v] = PETSC_TRUE;
573: thisIdx.push_back(g.local_to_global(v));
574: thisSize++;
575: for (i=cone_succ[v].begin();i!=cone_succ[v].end();i++) {
576: q.push(*i);
577: }
578: for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
579: PetscInt t = target(*e,g);
580: if (!taken[t]) {
581: mycone_pq.push(PQNode(t,v,get(edge_length_g,*e)));
582: boundary += get(edge_weight_g,*e);
583: edgeCount++;
584: } else if (mycone_taken[t]) {
585: boundary -= get(edge_weight_g,*e);
586: initialInternalConeEdges++;
587: }
588: }
589: }
590: if (initialInternalConeEdges < edgesLeft) {
591: while (initialInternalConeEdges == 0 ?
592: boundary > (edgeCount+1)*log(edgesLeft+1)*2.0/(log(2.0)*epsilon*radius) :
593: boundary > (edgeCount)*log(edgesLeft*1.0/initialInternalConeEdges)*2.0/(log(2.0)*epsilon*radius))
594: {
595: assert(!mycone_pq.empty());
596: node = mycone_pq.top();mycone_pq.pop();
597: if (!mycone_taken[node.vertex]) {
598: q.push(node.vertex);
599: while (!q.empty()) {
600: v = q.front();q.pop();
601: taken[v] = PETSC_TRUE;
602: mycone_taken[v] = PETSC_TRUE;
603: thisIdx.push_back(g.local_to_global(v));
604: thisSize++;
605: for (i=cone_succ[v].begin();i!=cone_succ[v].end();i++) {
606: q.push(*i);
607: }
608: for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
609: PetscInt t = target(*e,g);
610: if (!taken[t]) {
611: mycone_pq.push(PQNode(t,v,node.dist+get(edge_length_g,*e)));
612: boundary += get(edge_weight_g,*e);
613: edgeCount++;
614: } else if (mycone_taken[t]) {
615: boundary -= get(edge_weight_g,*e);
616: }
617: }
618: }
619: }
620: }
621: }
622: edgesLeft -= edgeCount;
623: size.push_back(thisSize);
624: idx.push_back(thisIdx);
625: x.push_back(anchor);
626: y.push_back(pred[anchor]);
627: k++;
628: }
629: }
630:
631:
633: /*
634: // pseudo cone cut
635: while (!pq.empty()) {
636: node = pq.top();pq.pop();
638: PetscInt thisSize = 1;
639: std::vector<PetscInt> thisIdx;
640: std::queue<PetscInt> q;
642: thisIdx.push_back(g.local_to_global(node.vertex));
643: for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
644: q.push(*i);
645: }
647: PetscInt v;
648: while (!q.empty()) {
649: v = q.front();q.pop();
650: thisSize++;
651: thisIdx.push_back(g.local_to_global(v));
652: for (i=succ[v].begin();i!=succ[v].end();i++) {
653: q.push(*i);
654: }
655: }
656: size.push_back(thisSize);
657: idx.push_back(thisIdx);
658: x.push_back(node.vertex);
659: y.push_back(pred[node.vertex]);
660: k++;
661: }
662: */
664:
667: return(0);
668: }
670: /* -------------------------------------------------------------------------- */
671: /*
672: AugmentSpanningTree
674: Input Parameters:
675: . g - input graph with spanning tree defined by vertex_parent and vertex_children;
676: vertex_depth gives distance in tree from root
677: . maxCong - target upper bound on congestion
679: Output Parameter:
680: . g - edges in augmented spanning tree are marked with edge_keep = PETSC_TRUE
681: . maxCong - an actual upper bound on congestion
682: */
685: PetscErrorCode AugmentSpanningTree(Graph& g,const PetscInt root,PetscScalar& maxCong)
686: {
687: // const PetscInt n = num_vertices(g);
688: const PetscInt m = num_edges(g);
689: //PetscInt i;
691: //PetscInt *component; // maps each vertex to a vertex component
692: //std::vector<PetscScalar>
693: // maxCongestion; /* maps each edge component to an upper bound on the
694: // congestion through any of its edges */
696: const EdgeIndex edge_index_g = get(edge_index_t(),g);
700: std::vector<Component*> componentList;
701: std::vector<ComponentPair> edgeComponentMap(m);
703: DecomposeSpanningTree(g,root,maxCong,maxCong,
704: componentList,edgeComponentMap);
705: /*
706: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"COMPONENTS:\n");
707: for (int i=0; i<componentList.size(); i++) {
708: for (int j=0; j<componentList[i]->vertices.size(); j++) {
709: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"%d ",componentList[i]->vertices[j]);
710: }
711: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"EDGES: ");
713: graph_traits<Graph>::edge_iterator e, e_end;
714: for (tie(e,e_end)=edges(g); e!=e_end; e++) {
715: if (edgeComponentMap[get(edge_index_g,*e)].getIndex(componentList[i]) != -1) {
716: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"(%d,%d) ",source(*e,g),target(*e,g));
717: }
718: }
719: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"\n");
720: }
721: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"END COMPONENTS\n");
722: */
723: AddBridges(g,componentList,edgeComponentMap,maxCong);
725: return(0);
726: }
730: /* -------------------------------------------------------------------------- */
731: /*
732: DecomposeSpanningTree
734: Input Parameters:
735: . g - input graph with spanning tree defined by vertex_parent and vertex_children;
736: vertex_depth gives distance in tree from root
737: . component - a preallocated array of length num_vertices(g)
738: . edgeComponent - a preallocated vector of length num_edges(g)
740: Output Parameter:
741: . component - a vector mapping each vertex to a component
742: */
745: PetscErrorCode DecomposeSpanningTree(Graph& g,const PetscInt root,
746: const PetscScalar maxInternalStretch,
747: const PetscScalar maxExternalWeight,
748: std::vector<Component*>& componentList,
749: std::vector<ComponentPair>& edgeComponentMap)
750: {
755: Component* currComponent;
756: PetscScalar currInternalStretch, currExternalWeight;
758: //Component::next_id = 0;
759: //Component::max_id = num_edges(g);
760: DecomposeSubTree(g,root,
761: maxInternalStretch,maxExternalWeight,
762: componentList,edgeComponentMap,
763: currComponent,
764: currInternalStretch,
765: currExternalWeight);
766:
767: componentList.push_back(currComponent);
769: return(0);
770: }
775: PetscErrorCode DecomposeSubTree(Graph& g,const PetscInt root,
776: const PetscScalar maxInternalStretch,
777: const PetscScalar maxExternalWeight,
778: std::vector<Component*>& componentList,
779: std::vector<ComponentPair>& edgeComponentMap,
780: Component*& currComponent,
781: PetscScalar& currInternalStretch,
782: PetscScalar& currExternalWeight)
783: {
784: const EdgeWeight edge_weight_g = get(edge_weight_t(),g);
785: const EdgeIndex edge_index_g = get(edge_index_t(),g);
786: const EdgeKeep edge_keep_g = get(edge_keep_t(),g);
787: const VertexParent vertex_parent_g = get(vertex_parent_t(),g);
788: const VertexChildren vertex_children_g = get(vertex_children_t(),g);
789: const VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
790: const PetscScalar rootDepth = get(vertex_depth_g,root);
791: std::vector<PetscInt>::const_iterator i,j;
792: graph_traits<Graph>::out_edge_iterator e, e_end;
794: PetscScalar newInternalStretch, newExternalWeight;
795: PetscInt v,w,edgeIndex,compIndex;
796: PetscScalar weight;
800: currComponent = new Component();
801: currComponent->vertices.push_back(root);
802: currInternalStretch = 0;
803: currExternalWeight = 0;
804: // temporarily add all external edges to component, but may remove some at end
805: for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
806: if (!get(edge_keep_g,*e)) {
807: edgeIndex = get(edge_index_g,*e);
808:
809: if (edgeComponentMap[edgeIndex].get(0) == PETSC_NULL) {
810: edgeComponentMap[edgeIndex].put(0,currComponent);
811: } else {
812: assert(edgeComponentMap[edgeIndex].get(1) == PETSC_NULL);
813: edgeComponentMap[edgeIndex].put(1,currComponent);
814: }
815: }
816: }
818: std::vector<PetscInt> children = get(vertex_children_g,root);
820: Component *childComponent;
821: PetscScalar childInternalStretch, childExternalWeight;
822: for (i=children.begin();i!=children.end();i++) {
823: PetscInt child = *i;
824: DecomposeSubTree(g,child,maxInternalStretch,maxExternalWeight,
825: componentList,edgeComponentMap,
826: childComponent,
827: childInternalStretch,childExternalWeight);
829: newInternalStretch = currInternalStretch + childInternalStretch;
830: newExternalWeight = currExternalWeight;
832: for (j = childComponent->vertices.begin(), v = *j;
833: j != childComponent->vertices.end() && (newInternalStretch <= maxInternalStretch);
834: v = *(++j)) {
835: for (tie(e,e_end)=out_edges(v,g);
836: e!=e_end && (newInternalStretch <= maxInternalStretch);
837: e++) {
838: if (!get(edge_keep_g,*e)) {
839: w = target(*e,g);
840: edgeIndex = get(edge_index_g,*e);
841: compIndex = edgeComponentMap[edgeIndex].getIndex(childComponent);
842:
843: if (compIndex != -1) {
844: weight = get(edge_weight_g,*e);
845:
846: if (edgeComponentMap[edgeIndex].get(1-compIndex) == currComponent) {
847: newExternalWeight -= weight;
848: newInternalStretch +=
849: (get(vertex_depth_g,v) + get(vertex_depth_g,w) - 2*rootDepth) * weight;
850: } else {
851: newExternalWeight += weight;
852: }
853: }
854: }
855: }
856: }
858: if (newInternalStretch <= maxInternalStretch && newExternalWeight <= maxExternalWeight) {
859: // merge the components
861: currInternalStretch = newInternalStretch;
862: currExternalWeight = newExternalWeight;
864: for (j = childComponent->vertices.begin(), v = *j;
865: j != childComponent->vertices.end();
866: v = *(++j)) {
867: currComponent->vertices.push_back(v);
868: for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
869: if (!get(edge_keep_g,*e)) {
870: edgeIndex = get(edge_index_g,*e);
871: if (edgeComponentMap[edgeIndex].get(0) == childComponent) {
872: edgeComponentMap[edgeIndex].put(0,currComponent);
873: }
874: if (edgeComponentMap[edgeIndex].get(1) == childComponent) {
875: edgeComponentMap[edgeIndex].put(1,currComponent);
876: }
877: }
878: }
879: }
880: delete childComponent;
881: } else {
882: componentList.push_back(childComponent);
883: }
884: }
886: const Component *origCurrComponent = currComponent;
887: for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
888: edgeIndex = get(edge_index_g,*e);
889: if (!get(edge_keep_g,*e)) {
890: if (edgeComponentMap[edgeIndex].get(0) == origCurrComponent) {
891: compIndex = 0;
892: } else {
893: assert(edgeComponentMap[edgeIndex].get(1) == origCurrComponent);
894: compIndex = 1;
895: }
896:
897: if (edgeComponentMap[edgeIndex].get(1-compIndex) != origCurrComponent) {
898: weight = get(edge_weight_g,*e);
899: if (currExternalWeight + weight <= maxExternalWeight) {
900: currExternalWeight += weight;
901: } else {
902: componentList.push_back(currComponent);
903: currComponent = new Component();
904: currComponent->vertices.push_back(root);
905: currInternalStretch = 0;
906: currExternalWeight = 0;
907: }
908: edgeComponentMap[edgeIndex].put(compIndex,currComponent);
909: }
910: }
911: }
913: return(0);
914: }
919: PetscErrorCode AddBridges(Graph& g,
920: std::vector<Component*>& componentList,
921: std::vector<ComponentPair>& edgeComponentMap,
922: PetscScalar& maxCong) {
924: const PetscInt m = num_edges(g);
925: std::vector<PetscTruth> edgeSupported(m); // edgeSupported[i] if edge i's component pair has been bridged
926: const EdgeLength edge_length_g = get(edge_length_t(),g);
927: const EdgeWeight edge_weight_g = get(edge_weight_t(),g);
928: const EdgeIndex edge_index_g = get(edge_index_t(),g);
929: const EdgeKeep edge_keep_g = get(edge_keep_t(),g);
930: const VertexParent vertex_parent_g = get(vertex_parent_t(),g);
931: const VertexChildren vertex_children_g = get(vertex_children_t(),g);
932: const VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
933: PetscInt edgeIndex, eeIndex;
934: Component *comp1, *comp2;
935: PetscInt comp1size, comp2size, i, v, w, parent;
936: graph_traits<Graph>::edge_iterator e, e_end;
937: graph_traits<Graph>::out_edge_iterator ee, ee_end;
938: PetscScalar realMaxCong;
942: realMaxCong = 0;
944: for (tie(e,e_end)=edges(g); e!=e_end; e++) {
945: if (!get(edge_keep_g,*e)) {
946: edgeIndex = get(edge_index_g,*e);
947: comp1 = edgeComponentMap[edgeIndex].get(0);
948: comp2 = edgeComponentMap[edgeIndex].get(1);
949: if ((comp1 != comp2) && !edgeSupported[edgeIndex]) {
950: comp1size = comp1->vertices.size();
951: comp2size = comp2->vertices.size();
952: std::map<PetscInt,PetscScalar> congestionBelow1,weightBelow1;
953: std::map<PetscInt,PetscScalar> congestionBelow2,weightBelow2;
954: for (i=0; i<comp1size; i++) {
955: congestionBelow1[comp1->vertices[i]] = 0;
956: weightBelow1[comp1->vertices[i]] = 0;
957: }
958: for (i=0; i<comp2size; i++) {
959: congestionBelow2[comp2->vertices[i]] = 0;
960: weightBelow2[comp2->vertices[i]] = 0;
961: }
963: for (i=comp1size-1; i>=0; i--) {
964: v = comp1->vertices[i];
965: for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
966: if (!get(edge_keep_g,*ee)) {
967: eeIndex = get(edge_index_g,*ee);
968: if (edgeComponentMap[eeIndex].match(comp1,comp2) &&
969: weightBelow2.count(target(*ee,g)) > 0) {
970: edgeSupported[eeIndex] = PETSC_TRUE;
971: weightBelow1[v] += get(edge_weight_g,*ee);
972: }
973: }
974: }
975: if (i>0) {
976: parent = get(vertex_parent_g,v);
977: weightBelow1[parent] += weightBelow1[v];
978: congestionBelow1[parent] +=
979: weightBelow1[v]*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
980: }
981: }
982: for (i=1; i<comp1size; i++) {
983: v = comp1->vertices[i];
984: parent = get(vertex_parent_g,v);
985: congestionBelow1[v] = congestionBelow1[parent] -
986: (weightBelow1[comp1->vertices[0]] - 2*weightBelow1[v])*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
987: }
988:
989: for (i=comp2size-1; i>=0; i--) {
990: v = comp2->vertices[i];
991: for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
992: if (!get(edge_keep_g,*ee)) {
993: eeIndex = get(edge_index_g,*ee);
994: if (edgeComponentMap[eeIndex].match(comp1,comp2) &&
995: weightBelow1.count(target(*ee,g)) > 0) {
996: assert(edgeSupported[eeIndex] == PETSC_TRUE);
997: weightBelow2[v] += get(edge_weight_g,*ee);
998: }
999: }
1000: }
1001: if (i>0) {
1002: parent = get(vertex_parent_g,v);
1003: weightBelow2[parent] += weightBelow2[v];
1004: congestionBelow2[parent] +=
1005: weightBelow2[v]*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
1006: }
1007: }
1008: for (i=1; i<comp2size; i++) {
1009: v = comp2->vertices[i];
1010: parent = get(vertex_parent_g,v);
1011: congestionBelow2[v] = congestionBelow2[parent] -
1012: (weightBelow2[comp2->vertices[0]] - 2*weightBelow2[v])*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
1013: }
1014: /*
1015: for (std::map<PetscInt,PetscScalar>::iterator it = congestionBelow1.begin();
1016: it != congestionBelow1.end();
1017: it++) {
1018: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"congestionBelow1[%d]=%f\n",(*it).first,(*it).second);
1019: }
1020: for (std::map<PetscInt,PetscScalar>::iterator it = weightBelow1.begin();
1021: it != weightBelow1.end();
1022: it++) {
1023: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"weightBelow1[%d]=%f\n",(*it).first,(*it).second);
1024: }
1025: for (std::map<PetscInt,PetscScalar>::iterator it = congestionBelow2.begin();
1026: it != congestionBelow2.end();
1027: it++) {
1028: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"congestionBelow2[%d]=%f\n",(*it).first,(*it).second);
1029: }
1030: for (std::map<PetscInt,PetscScalar>::iterator it = weightBelow2.begin();
1031: it != weightBelow2.end();
1032: it++) {
1033: PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"weightBelow2[%d]=%f\n",(*it).first,(*it).second);
1034: }
1035: */
1036: assert(weightBelow1[comp1->vertices[0]] - weightBelow2[comp2->vertices[0]] < 1e-9 &&
1037: weightBelow1[comp1->vertices[0]] - weightBelow2[comp2->vertices[0]] > -1e-9);
1039: Edge bestEdge;
1040: PetscScalar bestCongestion = -1;
1041: for (i=0; i<comp1size; i++) {
1042: v = comp1->vertices[i];
1043: for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
1044: if (!get(edge_keep_g,*ee)) {
1045: eeIndex = get(edge_index_g,*ee);
1046: if (edgeComponentMap[eeIndex].match(comp1,comp2)) {
1047: w = target(*ee,g);
1048: PetscScalar newCongestion =
1049: weightBelow1[comp1->vertices[0]] * get(edge_length_g,*ee) +
1050: congestionBelow1[v] + congestionBelow2[w];
1051: if (bestCongestion < 0 || newCongestion < bestCongestion) {
1052: bestEdge = *ee;
1053: bestCongestion = newCongestion;
1054: }
1055: }
1056: }
1057: }
1058: }
1059: put(edge_keep_g,bestEdge,PETSC_TRUE);
1060: if (bestCongestion > realMaxCong)
1061: realMaxCong = bestCongestion;
1062: }
1063: }
1064: }
1065: maxCong = realMaxCong;
1067: return(0);
1068: }