Actual source code: lowstretch.cpp

  1: #define PETSCKSP_DLL

  3: #include <math.h>
  4: #include <queue>
  5: #include <private/pcimpl.h>   /*I "petscpc.h" I*/
  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, PetscBool,
 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:   PetscReal dist;

115:   PQNode() {}

117:   PQNode(const PetscInt v,const PetscReal d) {
118:     vertex = v;
119:     dist = d;
120:   }

122:   PQNode(const PetscInt v,const PetscInt p,const PetscReal 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:                                       PetscReal tol,PetscReal& maxCong);
139: PetscErrorCode AugmentedLowStretchSpanningTree(Mat mat,Mat *pre,PetscBool augment,
140:                                                PetscReal tol,PetscReal& 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:                                       PetscReal tol,PetscReal& 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,PetscBool augment,
207:                                                PetscReal tol,PetscReal& maxCong)
208: {
209: #ifndef PETSC_USE_COMPLEX
210:   PetscInt          *idx;
211:   PetscInt          start,end,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;
221: #endif
222:   PetscInt          n;
223:   PetscErrorCode    ierr;

226:   MatGetSize(mat,PETSC_NULL,&n);

228:   Graph g(n);

230:   EdgeKeep edge_keep_g = get(edge_keep_t(),g);

232: #ifdef PETSC_USE_COMPLEX
233:   SETERRQ(PETSC_ERR_SUP, "Complex numbers not supported for support graph PC");
234: #else
235:   PetscMalloc3(n,PetscScalar,&diag,n,PetscInt,&dnz,n,PetscInt,&onz);
236:   PetscMemzero(dnz, n * sizeof(PetscInt));
237:   PetscMemzero(onz, n * sizeof(PetscInt));
238:   MatGetOwnershipRange(mat, &start, &end);
239:   for (i=0; i<n; i++) {
240:     MatGetRow(mat,i,&ncols,&cols_c,&vals_c);
241:     diag[i] = 0;
242:     for (k=0; k<ncols; k++) {
243:       if (cols_c[k] == i) {
244:         diag[i] += vals_c[k];
245:       } else if (vals_c[k] != 0) {
246:         absval = vals_c[k]>0?vals_c[k]:-vals_c[k];
247:         diag[i] -= absval;
248:         if (cols_c[k] > i && absval > tol) {
249:           // we set edge_weight = absval; edge_length = 1.0/absval; edge_keep = PETSC_FALSE
250:           add_edge(i,cols_c[k],EDGE_PROPERTY(absval,1.0/absval,PETSC_FALSE),g);
251:         }
252:       }
253:     }
254:     MatRestoreRow(mat,i,&ncols,&cols_c,&vals_c);
255:   }

257:   PetscMalloc(n*sizeof(PetscInt),&idx);
258:   put(get(vertex_depth_t(),g),root,0);
259:   LowStretchSpanningTreeHelper(g,root,log(4.0/3.0)/(2.0*log((double)n)),idx);

261:   if (augment) {
262:     AugmentSpanningTree(g,root,maxCong);
263:   }

265:   pre = prefact;
266:   MatCreate(PETSC_COMM_WORLD,pre);
267:   MatSetSizes(*pre,PETSC_DECIDE,PETSC_DECIDE,n,n);
268:   MatSetType(*pre,MATAIJ);

270:   PetscMalloc3(1,PetscInt,&rows,n,PetscInt,&cols,n,PetscScalar,&vals);
271:   for (i=0; i<n; i++) {
272:     for (tie(e, e_end) = out_edges(i,g); e != e_end; e++) {
273:       if (get(edge_keep_g,*e)) {
274:         const PetscInt col =  target(*e,g);

276:         if (col >= start && col < end) {
277:           dnz[i]++;
278:         } else {
279:           onz[i]++;
280:         }
281:       }
282:     }
283:   }
284:   MatSeqAIJSetPreallocation(*pre, 0, dnz);
285:   MatMPIAIJSetPreallocation(*pre, 0, dnz, 0, onz);
286:   for (i=0; i<n; i++) {
287:     rows[0] = i;
288:     k = 0;
289:     for (tie(e, e_end) = out_edges(i,g); e != e_end; e++) {
290:       if (get(edge_keep_g,*e)) {
291:         cols[k++] = target(*e,g);
292:       }
293:     }
294:     MatGetValues(mat,1,rows,k,cols,vals);
295:     for (j=0; j<k; j++) {
296:       absval = vals[j]>0?vals[j]:-vals[j];
297:       diag[i] += absval;
298:     }
299:     cols[k] = i;
300:     vals[k] = diag[i];
301:     MatSetValues(*pre,1,rows,k+1,cols,vals,INSERT_VALUES);
302:   }
303:   PetscFree3(rows,cols,vals);
304:   PetscFree3(diag,dnz,onz);

306:   MatAssemblyBegin(*pre,MAT_FINAL_ASSEMBLY);
307:   MatAssemblyEnd(*pre,MAT_FINAL_ASSEMBLY);

309:   /*
310:     ISCreateGeneral(PETSC_COMM_WORLD,n,idx,&perm);
311:     ISSetPermutation(perm);
312:     ISInvertPermutation(perm,PETSC_DECIDE,&iperm);
313:     PetscFree(idx);
314:     ISView(perm,PETSC_VIEWER_STDOUT_SELF);
315:   */
316:   IS rperm, cperm;
317:   MatGetOrdering(*pre,MATORDERINGQMD,&rperm,&cperm);

319:   MatFactorInfoInitialize(&info);
320:   /*
321:   MatCholeskyFactorSymbolic(*pre,iperm,&info,prefact);
322:   MatCholeskyFactorNumeric(*pre,&info,prefact);
323:   MatDestroy(*pre);
324:   ISDestroy(perm);
325:   ISDestroy(iperm);
326:   */
327:   MatLUFactor(*pre,rperm,cperm,&info);
328:   ISDestroy(&rperm);
329:   ISDestroy(&cperm);
330: #endif
331:   return(0);
332: }

334: /* -------------------------------------------------------------------------- */
335: /*
336:    LowStretchSpanningTreeHelper

338:    Input Parameters:
339: .  g - input graph; all edges have edge_keep = PETSC_FALSE;
340:                     root has vertex_depth set to distance from global root
341: .  root - root vertex
342: .  alpha - parameter
343: .  perm - preallocated array of size num_vertices(g) in which to store vertex ordering

345:    Output Parameter:
346: .  g - edges in low-stretch spanning tree are marked with edge_keep = PETSC_TRUE;
347:        also vertex_parent and vertex_children are set (vertex_parent of global root is undefined)
348:        and vertex_depth is set to be distance from global root (where weight on edge is inverse distance)
349: .  perm - list of vertices in which a vertex precedes its parent (with vertices referred to by global index)
350:  */
353: PetscErrorCode LowStretchSpanningTreeHelper(Graph& g,const PetscInt root,const PetscScalar alpha,PetscInt perm[])
354: {
355:   PetscInt n,i,j,k;
356:   std::vector<PetscInt> size,x,y;
357:   std::vector<std::vector<PetscInt> > idx;


362:   EdgeLength edge_length_g = get(edge_length_t(),g);
363:   EdgeKeep edge_keep_g = get(edge_keep_t(),g);
364:   VertexParent vertex_parent_g = get(vertex_parent_t(),g);
365:   VertexChildren vertex_children_g = get(vertex_children_t(),g);
366:   VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
367:   n = num_vertices(g);

369:   if (n > 2) {
370:     StarDecomp(g,root,1.0/3,alpha,k,size,idx,x,y);
371:     j = n - size[0];
372:     Graph& g1 = g.create_subgraph(idx[0].begin(),idx[0].end());
373:     LowStretchSpanningTreeHelper(g1,g1.global_to_local(g.local_to_global(root)),alpha,perm+j);
374:     for (i=1;i<=k;i++) {
375:       Edge e = edge(x[i-1],y[i-1],g).first;
376:       put(edge_keep_g,e,PETSC_TRUE);
377:       put(vertex_parent_g,x[i-1],g.local_to_global(y[i-1]));
378:       get(vertex_children_g,y[i-1]).push_back(g.local_to_global(x[i-1]));
379:       put(vertex_depth_g,x[i-1],get(vertex_depth_g,y[i-1])+get(edge_length_g,e));

381:       j -= size[i];
382:       Graph& g1 = g.create_subgraph(idx[i].begin(),idx[i].end());
383:       LowStretchSpanningTreeHelper(g1,g1.global_to_local(g.local_to_global(x[i-1])),alpha,perm+j);
384:     }
385:   } else if (n == 2) {
386:     Edge e = *(out_edges(root,g).first);
387:     PetscInt t = target(e,g);

389:     put(edge_keep_g,e,PETSC_TRUE);
390:     put(vertex_parent_g,t,g.local_to_global(root));
391:     get(vertex_children_g,root).push_back(g.local_to_global(t));
392:     put(vertex_depth_g,t,get(vertex_depth_g,root)+get(edge_length_g,e));

394:     perm[0] = g.local_to_global(t);
395:     perm[1] = g.local_to_global(root);
396:   } else /* n == 1 */ {
397:     perm[0] = g.local_to_global(root);
398:   }

400:   return(0);
401: }



405: /* -------------------------------------------------------------------------- */
406: /*
407:    StarDecomp - calculate a star decomposition of the graph

409:    Input Parameters:
410: .  g - input graph
411: .  root - center of star decomposition
412: .  delta, epsilon - parameters of star decomposition

414:    Output Parameter:
415: .  k - number of partitions, not-including central one.
416: .  size[i] - number of vertices in partition i
417: .  idx[i] - list of vertices in partition i; partition 0 contains root
418: .  x[i-1] - "anchor" vertex of non-central partition i
419: .  y[i-i] - vertex in partition 0 forming "bridge" with x[i-1]
420:  */
423: PetscErrorCode StarDecomp(Graph g,const PetscInt root,const PetscScalar delta,const PetscScalar epsilon,
424:                           PetscInt& k,std::vector<PetscInt>& size,std::vector<std::vector<PetscInt> >& idx,
425:                           std::vector<PetscInt>& x,std::vector<PetscInt>& y)
426: {
427: #ifndef PETSC_USE_COMPLEX
428:   PetscInt n,m,edgesLeft;
430:   ShortestPathPriorityQueue pq;
431:   PetscScalar radius;
432:   PetscInt centerSize;
433:   std::vector<PetscInt> centerIdx;
434:   PQNode node;
435: #endif

438: #ifdef PETSC_USE_COMPLEX
439:   SETERRQ(PETSC_ERR_SUP, "Complex numbers not supported for support graph PC");
440: #else
441:   EdgeWeight edge_weight_g = get(edge_weight_t(),g);
442:   EdgeLength edge_length_g = get(edge_length_t(),g);
443:   n = num_vertices(g);
444:   m = num_edges(g);
445:   edgesLeft = m;

447:   std::vector<PetscInt> pred(n,-1);
448:   std::vector<PetscInt> *succ = new std::vector<PetscInt>[n];
449:   std::vector<PetscInt>::iterator i;
450:   PetscScalar *dist = new PetscScalar[n];
451:   std::vector<PetscBool> taken(n,PETSC_FALSE);

453:   /** form tree of shortest paths to root **/
454:   graph_traits<Graph>::out_edge_iterator e, e_end;
455:   for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
456:     PetscInt t = target(*e,g);
457:     pq.push(PQNode(t,root,get(edge_length_g,*e)));
458:   }
459:   pred[root] = root;
460:   while (!pq.empty()) {
461:     node = pq.top();pq.pop();
462:     if (pred[node.vertex] == -1) {
463:       succ[node.pred].push_back(node.vertex);
464:       pred[node.vertex] = node.pred;
465:       dist[node.vertex] = node.dist;
466:       for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
467:         PetscInt t = target(*e,g);
468:         if (pred[t] == -1) {
469:           pq.push(PQNode(t,node.vertex,node.dist+get(edge_length_g,*e)));
470:         }
471:       }
472:       radius = node.dist;
473:     }
474:   }

476:   /** BALL CUT **/
477:   for (i=succ[root].begin();i!=succ[root].end();i++) {
478:     pq.push(PQNode(*i,dist[*i]));
479:   }
480:   PetscScalar boundary = 0;
481:   PetscInt edgeCount = 0;
482:   centerIdx.push_back(g.local_to_global(root));
483:   taken[root] = PETSC_TRUE;
484:   centerSize = 1;
485:   for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
486:     boundary += get(edge_weight_g,*e);
487:     edgeCount++;
488:   }
489:   const PetscScalar minRadius = delta*radius;
490:   while (dist[pq.top().vertex] < minRadius) {
491:     assert(!pq.empty());
492:     node = pq.top();pq.pop();
493:     centerIdx.push_back(g.local_to_global(node.vertex));
494:     taken[node.vertex] = PETSC_TRUE;
495:     centerSize++;
496:     for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
497:       if (taken[target(*e,g)]) {
498:         boundary -= get(edge_weight_g,*e);
499:       } else {
500:         boundary += get(edge_weight_g,*e);
501:         edgeCount++;
502:       }
503:     }
504:     for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
505:       pq.push(PQNode(*i,dist[*i]));
506:     }
507:   }
508:   while (boundary > (edgeCount+1)*log((double)m)/(log(2.0)*(1.0-2.0*delta)*radius)) {
509:     assert(!pq.empty());
510:     node = pq.top();pq.pop();
511:     centerIdx.push_back(g.local_to_global(node.vertex));
512:     taken[node.vertex] = PETSC_TRUE;
513:     centerSize++;
514:     for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
515:       if (taken[target(*e,g)]) {
516:         boundary -= get(edge_weight_g,*e);
517:       } else {
518:         boundary += get(edge_weight_g,*e);
519:         edgeCount++;
520:       }
521:     }
522:     for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
523:       pq.push(PQNode(*i,dist[*i]));
524:     }
525:   }
526:   size.push_back(centerSize);
527:   idx.push_back(centerIdx);
528:   edgesLeft -= edgeCount;

530:   k = 0;
531:   assert(!pq.empty());
532:   std::queue<PetscInt> anchor_q;
533:   ShortestPathPriorityQueue cone_pq;
534:   std::vector<PetscInt> *cone_succ = new std::vector<PetscInt>[n];
535:   std::vector<PetscBool> cone_found(n,PETSC_FALSE);

537:   /** form tree of shortest paths to an anchor **/
538:   while (!pq.empty()) {
539:     node = pq.top();pq.pop();
540:     cone_found[node.vertex] = PETSC_TRUE;
541:     anchor_q.push(node.vertex);
542:     for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
543:       PetscInt t = target(*e,g);
544:       if (!taken[t]) {
545:         cone_pq.push(PQNode(t,node.vertex,get(edge_length_g,*e)));
546:       }
547:     }
548:   }
549:   while (!cone_pq.empty()) {
550:     node = cone_pq.top();cone_pq.pop();
551:     if (!cone_found[node.vertex]) {
552:       cone_succ[node.pred].push_back(node.vertex);
553:       cone_found[node.vertex] = PETSC_TRUE;
554:       for (tie(e,e_end)=out_edges(node.vertex,g); e!=e_end; e++) {
555:         PetscInt t = target(*e,g);
556:         if (!taken[t] && !cone_found[t]) {
557:           cone_pq.push(PQNode(t,node.vertex,node.dist+get(edge_length_g,*e)));
558:         }
559:       }
560:     }
561:   }

563:   while (!anchor_q.empty()) {
564:     /** CONE CUT **/
565:     PetscInt anchor = anchor_q.front();anchor_q.pop();
566:     if (!taken[anchor]) {
567:       PetscInt v;
568:       PetscInt thisSize = 0;
569:       std::vector<PetscInt> thisIdx;
570:       std::queue<PetscInt> q;
571:       ShortestPathPriorityQueue mycone_pq;
572:       std::vector<PetscBool> mycone_taken(n,PETSC_FALSE);
573:       PetscInt initialInternalConeEdges = 0;

575:       boundary = 0;
576:       edgeCount = 0;
577:       q.push(anchor);
578:       while (!q.empty()) {
579:         v = q.front();q.pop();
580:         taken[v] = PETSC_TRUE;
581:         mycone_taken[v] = PETSC_TRUE;
582:         thisIdx.push_back(g.local_to_global(v));
583:         thisSize++;
584:         for (i=cone_succ[v].begin();i!=cone_succ[v].end();i++) {
585:           q.push(*i);
586:         }
587:         for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
588:           PetscInt t = target(*e,g);
589:           if (!taken[t]) {
590:             mycone_pq.push(PQNode(t,v,get(edge_length_g,*e)));
591:             boundary += get(edge_weight_g,*e);
592:             edgeCount++;
593:           } else if (mycone_taken[t]) {
594:             boundary -= get(edge_weight_g,*e);
595:             initialInternalConeEdges++;
596:           }
597:         }
598:       }
599:       if (initialInternalConeEdges < edgesLeft) {
600:         while (initialInternalConeEdges == 0 ?
601:                boundary > (edgeCount+1)*log((double)(edgesLeft+1))*2.0/(log(2.0)*epsilon*radius) :
602:                boundary > (edgeCount)*log(edgesLeft*1.0/initialInternalConeEdges)*2.0/(log(2.0)*epsilon*radius))
603:           {
604:             assert(!mycone_pq.empty());
605:             node = mycone_pq.top();mycone_pq.pop();
606:             if (!mycone_taken[node.vertex]) {
607:               q.push(node.vertex);
608:               while (!q.empty()) {
609:                 v = q.front();q.pop();
610:                 taken[v] = PETSC_TRUE;
611:                 mycone_taken[v] = PETSC_TRUE;
612:                 thisIdx.push_back(g.local_to_global(v));
613:                 thisSize++;
614:                 for (i=cone_succ[v].begin();i!=cone_succ[v].end();i++) {
615:                   q.push(*i);
616:                 }
617:                 for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
618:                   PetscInt t = target(*e,g);
619:                   if (!taken[t]) {
620:                     mycone_pq.push(PQNode(t,v,node.dist+get(edge_length_g,*e)));
621:                     boundary += get(edge_weight_g,*e);
622:                     edgeCount++;
623:                   } else if (mycone_taken[t]) {
624:                     boundary -= get(edge_weight_g,*e);
625:                   }
626:                 }
627:               }
628:             }
629:           }
630:       }
631:       edgesLeft -= edgeCount;
632:       size.push_back(thisSize);
633:       idx.push_back(thisIdx);
634:       x.push_back(anchor);
635:       y.push_back(pred[anchor]);
636:       k++;
637:     }
638:   }
639: #endif    
640: 

642:   /*
643:   // pseudo cone cut
644:   while (!pq.empty()) {
645:     node = pq.top();pq.pop();

647:     PetscInt thisSize = 1;
648:     std::vector<PetscInt> thisIdx;
649:     std::queue<PetscInt> q;

651:     thisIdx.push_back(g.local_to_global(node.vertex));
652:     for (i=succ[node.vertex].begin();i!=succ[node.vertex].end();i++) {
653:       q.push(*i);
654:     }

656:     PetscInt v;
657:     while (!q.empty()) {
658:       v = q.front();q.pop();
659:       thisSize++;
660:       thisIdx.push_back(g.local_to_global(v));
661:       for (i=succ[v].begin();i!=succ[v].end();i++) {
662:         q.push(*i);
663:       }
664:     }
665:     size.push_back(thisSize);
666:     idx.push_back(thisIdx);
667:     x.push_back(node.vertex);
668:     y.push_back(pred[node.vertex]);
669:     k++;
670:   }
671:   */
672: #ifndef PETSC_USE_COMPLEX
673:   delete [] succ;
674:   delete [] dist;
675:   delete [] cone_succ;
676: #endif
677:   return(0);
678: }

680: /* -------------------------------------------------------------------------- */
681: /*
682:    AugmentSpanningTree

684:    Input Parameters:
685: .  g - input graph with spanning tree defined by vertex_parent and vertex_children;
686:        vertex_depth gives distance in tree from root
687: .  maxCong - target upper bound on congestion

689:    Output Parameter:
690: .  g - edges in augmented spanning tree are marked with edge_keep = PETSC_TRUE
691: .  maxCong - an actual upper bound on congestion
692:  */
695: PetscErrorCode AugmentSpanningTree(Graph& g,const PetscInt root,PetscScalar& maxCong)
696: {
697:   //  const PetscInt n = num_vertices(g);
698:   const PetscInt m = num_edges(g);
699:   //PetscInt i;
701:   //PetscInt *component;  // maps each vertex to a vertex component
702:   //std::vector<PetscScalar>
703:   //  maxCongestion;       /* maps each edge component to an upper bound on the
704:   //                            congestion through any of its edges */

706:   const EdgeIndex edge_index_g = get(edge_index_t(),g);


710:   std::vector<Component*> componentList;
711:   std::vector<ComponentPair> edgeComponentMap(m);

713:   DecomposeSpanningTree(g,root,maxCong,maxCong,
714:                                componentList,edgeComponentMap);
715:   /*
716:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"COMPONENTS:\n");
717:   for (int i=0; i<componentList.size(); i++) {
718:     for (int j=0; j<componentList[i]->vertices.size(); j++) {
719:       PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"%d ",componentList[i]->vertices[j]);
720:     }
721:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"EDGES: ");

723:     graph_traits<Graph>::edge_iterator e, e_end;  
724:     for (tie(e,e_end)=edges(g); e!=e_end; e++) {
725:       if (edgeComponentMap[get(edge_index_g,*e)].getIndex(componentList[i]) != -1) {
726:         PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"(%d,%d) ",source(*e,g),target(*e,g));
727:       }
728:     }
729:     PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"\n");
730:   }
731:   PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"END COMPONENTS\n");
732:   */
733:   AddBridges(g,componentList,edgeComponentMap,maxCong);

735:   return(0);
736: }



740: /* -------------------------------------------------------------------------- */
741: /*
742:    DecomposeSpanningTree

744:    Input Parameters:
745: .  g - input graph with spanning tree defined by vertex_parent and vertex_children;
746:        vertex_depth gives distance in tree from root
747: .  component - a preallocated array of length num_vertices(g)
748: .  edgeComponent - a preallocated vector of length num_edges(g)

750:    Output Parameter:
751: .  component - a vector mapping each vertex to a component
752:  */
755: PetscErrorCode DecomposeSpanningTree(Graph& g,const PetscInt root,
756:                                      const PetscScalar maxInternalStretch,
757:                                      const PetscScalar maxExternalWeight,
758:                                      std::vector<Component*>& componentList,
759:                                      std::vector<ComponentPair>& edgeComponentMap)
760: {


765:   Component* currComponent;
766:   PetscScalar currInternalStretch, currExternalWeight;

768:   //Component::next_id = 0;
769:   //Component::max_id = num_edges(g);
770:   DecomposeSubTree(g,root,
771:                           maxInternalStretch,maxExternalWeight,
772:                           componentList,edgeComponentMap,
773:                           currComponent,
774:                           currInternalStretch,
775:                           currExternalWeight);
776: 
777:   componentList.push_back(currComponent);

779:   return(0);
780: }


785: PetscErrorCode DecomposeSubTree(Graph& g,const PetscInt root,
786:                                 const PetscScalar maxInternalStretch,
787:                                 const PetscScalar maxExternalWeight,
788:                                 std::vector<Component*>& componentList,
789:                                 std::vector<ComponentPair>& edgeComponentMap,
790:                                 Component*& currComponent,
791:                                 PetscScalar& currInternalStretch,
792:                                 PetscScalar& currExternalWeight)
793: {
794: #ifndef PETSC_USE_COMPLEX
795:   const EdgeWeight edge_weight_g = get(edge_weight_t(),g);
796:   const EdgeIndex edge_index_g = get(edge_index_t(),g);
797:   const EdgeKeep edge_keep_g = get(edge_keep_t(),g);
798:   const VertexParent vertex_parent_g = get(vertex_parent_t(),g);
799:   const VertexChildren vertex_children_g = get(vertex_children_t(),g);
800:   const VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
801:   const PetscScalar rootDepth = get(vertex_depth_g,root);
802:   std::vector<PetscInt>::const_iterator i,j;
803:   graph_traits<Graph>::out_edge_iterator e, e_end;
805:   PetscScalar newInternalStretch, newExternalWeight;
806:   PetscInt v,w,edgeIndex,compIndex;
807:   PetscScalar weight;
808: #endif

811: #ifdef PETSC_USE_COMPLEX
812:   SETERRQ(PETSC_ERR_SUP, "Complex numbers not supported for support graph PC");
813: #else
814:   currComponent = new Component();
815:   currComponent->vertices.push_back(root);
816:   currInternalStretch = 0;
817:   currExternalWeight = 0;
818:   // temporarily add all external edges to component, but may remove some at end
819:   for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
820:     if (!get(edge_keep_g,*e)) {
821:       edgeIndex = get(edge_index_g,*e);
822: 
823:       if (edgeComponentMap[edgeIndex].get(0) == PETSC_NULL) {
824:         edgeComponentMap[edgeIndex].put(0,currComponent);
825:       } else {
826:         assert(edgeComponentMap[edgeIndex].get(1) == PETSC_NULL);
827:         edgeComponentMap[edgeIndex].put(1,currComponent);
828:       }
829:     }
830:   }

832:   std::vector<PetscInt> children = get(vertex_children_g,root);

834:   Component *childComponent;
835:   PetscScalar childInternalStretch, childExternalWeight;
836:   for (i=children.begin();i!=children.end();i++) {
837:     PetscInt child = *i;
838:     DecomposeSubTree(g,child,maxInternalStretch,maxExternalWeight,
839:                             componentList,edgeComponentMap,
840:                             childComponent,
841:                             childInternalStretch,childExternalWeight);

843:     newInternalStretch = currInternalStretch + childInternalStretch;
844:     newExternalWeight = currExternalWeight;

846:     for (j = childComponent->vertices.begin(), v = *j;
847:          j != childComponent->vertices.end() && (newInternalStretch <= maxInternalStretch);
848:          v = *(++j)) {
849:       for (tie(e,e_end)=out_edges(v,g);
850:            e!=e_end && (newInternalStretch <= maxInternalStretch);
851:            e++) {
852:         if (!get(edge_keep_g,*e)) {
853:           w = target(*e,g);
854:           edgeIndex = get(edge_index_g,*e);
855:           compIndex = edgeComponentMap[edgeIndex].getIndex(childComponent);
856: 
857:           if (compIndex != -1) {
858:             weight = get(edge_weight_g,*e);
859: 
860:             if (edgeComponentMap[edgeIndex].get(1-compIndex) == currComponent) {
861:               newExternalWeight -= weight;
862:               newInternalStretch +=
863:                 (get(vertex_depth_g,v) + get(vertex_depth_g,w) - 2*rootDepth) * weight;
864:             } else {
865:               newExternalWeight += weight;
866:             }
867:           }
868:         }
869:       }
870:     }

872:     if (newInternalStretch <= maxInternalStretch && newExternalWeight <= maxExternalWeight) {
873:       // merge the components

875:       currInternalStretch = newInternalStretch;
876:       currExternalWeight = newExternalWeight;

878:       for (j = childComponent->vertices.begin(), v = *j;
879:            j != childComponent->vertices.end();
880:            v = *(++j)) {
881:         currComponent->vertices.push_back(v);
882:         for (tie(e,e_end)=out_edges(v,g); e!=e_end; e++) {
883:           if (!get(edge_keep_g,*e)) {
884:             edgeIndex = get(edge_index_g,*e);
885:             if (edgeComponentMap[edgeIndex].get(0) == childComponent) {
886:               edgeComponentMap[edgeIndex].put(0,currComponent);
887:             }
888:             if (edgeComponentMap[edgeIndex].get(1) == childComponent) {
889:               edgeComponentMap[edgeIndex].put(1,currComponent);
890:             }
891:           }
892:         }
893:       }
894:       delete childComponent;
895:     } else {
896:       componentList.push_back(childComponent);
897:     }
898:   }

900:   const Component *origCurrComponent = currComponent;
901:   for (tie(e,e_end)=out_edges(root,g); e!=e_end; e++) {
902:     edgeIndex = get(edge_index_g,*e);
903:     if (!get(edge_keep_g,*e)) {
904:       if (edgeComponentMap[edgeIndex].get(0) == origCurrComponent) {
905:         compIndex = 0;
906:       } else {
907:         assert(edgeComponentMap[edgeIndex].get(1) == origCurrComponent);
908:         compIndex = 1;
909:       }
910: 
911:       if (edgeComponentMap[edgeIndex].get(1-compIndex) != origCurrComponent) {
912:         weight = get(edge_weight_g,*e);
913:         if (currExternalWeight + weight <= maxExternalWeight) {
914:           currExternalWeight += weight;
915:         } else {
916:           componentList.push_back(currComponent);
917:           currComponent = new Component();
918:           currComponent->vertices.push_back(root);
919:           currInternalStretch = 0;
920:           currExternalWeight = 0;
921:         }
922:         edgeComponentMap[edgeIndex].put(compIndex,currComponent);
923:       }
924:     }
925:   }
926: #endif
927:   return(0);
928: }


933: PetscErrorCode AddBridges(Graph& g,
934:                           std::vector<Component*>& componentList,
935:                           std::vector<ComponentPair>& edgeComponentMap,
936:                           PetscScalar& maxCong) {

938: #ifndef PETSC_USE_COMPLEX
939:   const PetscInt m = num_edges(g);
940:   std::vector<PetscBool> edgeSupported(m); // edgeSupported[i] if edge i's component pair has been bridged
941:   const EdgeLength edge_length_g = get(edge_length_t(),g);
942:   const EdgeWeight edge_weight_g = get(edge_weight_t(),g);
943:   const EdgeIndex edge_index_g = get(edge_index_t(),g);
944:   const EdgeKeep edge_keep_g = get(edge_keep_t(),g);
945:   const VertexParent vertex_parent_g = get(vertex_parent_t(),g);
946:   const VertexChildren vertex_children_g = get(vertex_children_t(),g);
947:   const VertexDepth vertex_depth_g = get(vertex_depth_t(),g);
948:   PetscInt edgeIndex, eeIndex;
949:   Component *comp1, *comp2;
950:   PetscInt comp1size, comp2size, i, v, w, parent;
951:   graph_traits<Graph>::edge_iterator e, e_end;
952:   graph_traits<Graph>::out_edge_iterator ee, ee_end;
953:   PetscScalar realMaxCong;
954: #endif

957: #ifdef PETSC_USE_COMPLEX
958:   SETERRQ(PETSC_ERR_SUP, "Complex numbers not supported for support graph PC");
959: #else
960:   realMaxCong = 0;

962:   for (tie(e,e_end)=edges(g); e!=e_end; e++) {
963:     if (!get(edge_keep_g,*e)) {
964:       edgeIndex = get(edge_index_g,*e);
965:       comp1 = edgeComponentMap[edgeIndex].get(0);
966:       comp2 = edgeComponentMap[edgeIndex].get(1);
967:       if ((comp1 != comp2) && !edgeSupported[edgeIndex]) {
968:         comp1size = comp1->vertices.size();
969:         comp2size = comp2->vertices.size();
970:         std::map<PetscInt,PetscScalar> congestionBelow1,weightBelow1;
971:         std::map<PetscInt,PetscScalar> congestionBelow2,weightBelow2;
972:         for (i=0; i<comp1size; i++) {
973:           congestionBelow1[comp1->vertices[i]] = 0;
974:           weightBelow1[comp1->vertices[i]] = 0;
975:         }
976:         for (i=0; i<comp2size; i++) {
977:           congestionBelow2[comp2->vertices[i]] = 0;
978:           weightBelow2[comp2->vertices[i]] = 0;
979:         }

981:         for (i=comp1size-1; i>=0; i--) {
982:           v = comp1->vertices[i];
983:           for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
984:             if (!get(edge_keep_g,*ee)) {
985:               eeIndex = get(edge_index_g,*ee);
986:               if (edgeComponentMap[eeIndex].match(comp1,comp2) &&
987:                   weightBelow2.count(target(*ee,g)) > 0) {
988:                 edgeSupported[eeIndex] = PETSC_TRUE;
989:                 weightBelow1[v] += get(edge_weight_g,*ee);
990:               }
991:             }
992:           }
993:           if (i>0) {
994:             parent = get(vertex_parent_g,v);
995:             weightBelow1[parent] += weightBelow1[v];
996:             congestionBelow1[parent] +=
997:               weightBelow1[v]*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
998:           }
999:         }
1000:         for (i=1; i<comp1size; i++) {
1001:           v = comp1->vertices[i];
1002:           parent = get(vertex_parent_g,v);
1003:           congestionBelow1[v] = congestionBelow1[parent] -
1004:             (weightBelow1[comp1->vertices[0]] - 2*weightBelow1[v])*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
1005:         }
1006: 
1007:         for (i=comp2size-1; i>=0; i--) {
1008:           v = comp2->vertices[i];
1009:           for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
1010:             if (!get(edge_keep_g,*ee)) {
1011:               eeIndex = get(edge_index_g,*ee);
1012:               if (edgeComponentMap[eeIndex].match(comp1,comp2) &&
1013:                   weightBelow1.count(target(*ee,g)) > 0) {
1014:                 assert(edgeSupported[eeIndex] == PETSC_TRUE);
1015:                 weightBelow2[v] += get(edge_weight_g,*ee);
1016:               }
1017:             }
1018:           }
1019:           if (i>0) {
1020:             parent = get(vertex_parent_g,v);
1021:             weightBelow2[parent] += weightBelow2[v];
1022:             congestionBelow2[parent] +=
1023:               weightBelow2[v]*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
1024:           }
1025:         }
1026:         for (i=1; i<comp2size; i++) {
1027:           v = comp2->vertices[i];
1028:           parent = get(vertex_parent_g,v);
1029:           congestionBelow2[v] = congestionBelow2[parent] -
1030:             (weightBelow2[comp2->vertices[0]] - 2*weightBelow2[v])*(get(vertex_depth_g,v)-get(vertex_depth_g,parent));
1031:         }
1032:         /*
1033:         for (std::map<PetscInt,PetscScalar>::iterator it = congestionBelow1.begin();
1034:              it != congestionBelow1.end();
1035:              it++) {
1036:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"congestionBelow1[%d]=%f\n",(*it).first,(*it).second);
1037:         }
1038:         for (std::map<PetscInt,PetscScalar>::iterator it = weightBelow1.begin();
1039:              it != weightBelow1.end();
1040:              it++) {
1041:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"weightBelow1[%d]=%f\n",(*it).first,(*it).second);
1042:         }
1043:         for (std::map<PetscInt,PetscScalar>::iterator it = congestionBelow2.begin();
1044:              it != congestionBelow2.end();
1045:              it++) {
1046:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"congestionBelow2[%d]=%f\n",(*it).first,(*it).second);
1047:         }
1048:         for (std::map<PetscInt,PetscScalar>::iterator it = weightBelow2.begin();
1049:              it != weightBelow2.end();
1050:              it++) {
1051:           PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_SELF,"weightBelow2[%d]=%f\n",(*it).first,(*it).second);
1052:         }
1053:         */
1054:         assert(weightBelow1[comp1->vertices[0]] - weightBelow2[comp2->vertices[0]] < 1e-9 &&
1055:                weightBelow1[comp1->vertices[0]] - weightBelow2[comp2->vertices[0]] > -1e-9);

1057:         Edge bestEdge;
1058:         PetscScalar bestCongestion = -1;
1059:         for (i=0; i<comp1size; i++) {
1060:           v = comp1->vertices[i];
1061:           for (tie(ee,ee_end)=out_edges(v,g); ee!=ee_end; ee++) {
1062:             if (!get(edge_keep_g,*ee)) {
1063:               eeIndex = get(edge_index_g,*ee);
1064:               if (edgeComponentMap[eeIndex].match(comp1,comp2)) {
1065:                 w = target(*ee,g);
1066:                 PetscScalar newCongestion =
1067:                   weightBelow1[comp1->vertices[0]] * get(edge_length_g,*ee) +
1068:                   congestionBelow1[v] + congestionBelow2[w];
1069:                 if (bestCongestion < 0 || newCongestion < bestCongestion) {
1070:                   bestEdge = *ee;
1071:                   bestCongestion = newCongestion;
1072:                 }
1073:               }
1074:             }
1075:           }
1076:         }
1077:         put(edge_keep_g,bestEdge,PETSC_TRUE);
1078:         if (bestCongestion > realMaxCong)
1079:           realMaxCong = bestCongestion;
1080:       }
1081:     }
1082:   }
1083:   maxCong = realMaxCong;
1084: #endif
1085:   return(0);
1086: }