1 /* -*- mode: C++; indent-tabs-mode: nil; -*-
3 * This file is a part of LEMON, a generic C++ optimization library.
5 * Copyright (C) 2003-2009
6 * Egervary Jeno Kombinatorikus Optimalizalasi Kutatocsoport
7 * (Egervary Research Group on Combinatorial Optimization, EGRES).
9 * Permission to use, modify and distribute this software is granted
10 * provided that this copyright notice appears in all copies. For
11 * precise terms see the accompanying LICENSE file.
13 * This software is provided "AS IS" with no warranty of any kind,
14 * express or implied, and with no claim as to its suitability for any
21 /// \brief Special plane digraph generator.
23 /// Graph generator application for various types of plane graphs.
29 /// for more info on the usage.
34 #include <lemon/list_graph.h>
35 #include <lemon/random.h>
36 #include <lemon/dim2.h>
37 #include <lemon/bfs.h>
38 #include <lemon/counter.h>
39 #include <lemon/suurballe.h>
40 #include <lemon/graph_to_eps.h>
41 #include <lemon/lgf_writer.h>
42 #include <lemon/arg_parser.h>
43 #include <lemon/euler.h>
44 #include <lemon/math.h>
45 #include <lemon/kruskal.h>
46 #include <lemon/time_measure.h>
48 using namespace lemon;
50 typedef dim2::Point<double> Point;
52 GRAPH_TYPEDEFS(ListGraph);
61 std::vector<Node> nodes;
62 ListGraph::NodeMap<Point> coords(g);
67 for(EdgeIt e(g);e!=INVALID;++e)
68 tlen+=sqrt((coords[g.v(e)]-coords[g.u(e)]).normSquare());
74 const double EPSILON=1e-8;
75 bool tsp_improve(Node u, Node v)
77 double luv=std::sqrt((coords[v]-coords[u]).normSquare());
82 for(IncEdgeIt e(g,v2);(n=g.runningNode(e))==u2;++e) { }
85 if(luv+std::sqrt((coords[v2]-coords[u2]).normSquare())-EPSILON>
86 std::sqrt((coords[u]-coords[u2]).normSquare())+
87 std::sqrt((coords[v]-coords[v2]).normSquare()))
89 g.erase(findEdge(g,u,v));
90 g.erase(findEdge(g,u2,v2));
100 bool tsp_improve(Node u)
102 for(IncEdgeIt e(g,u);e!=INVALID;++e)
103 if(tsp_improve(u,g.runningNode(e))) return true;
112 for(NodeIt n(g);n!=INVALID;++n)
113 if(tsp_improve(n)) b=true;
119 for(int i=0;i<N;i++) g.addEdge(nodes[i],nodes[(i+1)%N]);
128 Line(Point _a,Point _b) :a(_a),b(_b) {}
129 Line(Node _a,Node _b) : a(coords[_a]),b(coords[_b]) {}
130 Line(const Arc &e) : a(coords[g.source(e)]),b(coords[g.target(e)]) {}
131 Line(const Edge &e) : a(coords[g.u(e)]),b(coords[g.v(e)]) {}
134 inline std::ostream& operator<<(std::ostream &os, const Line &l)
136 os << l.a << "->" << l.b;
140 bool cross(Line a, Line b)
142 Point ao=rot90(a.b-a.a);
143 Point bo=rot90(b.b-b.a);
144 return (ao*(b.a-a.a))*(ao*(b.b-a.a))<0 &&
145 (bo*(a.a-b.a))*(bo*(a.b-b.a))<0;
155 bool pedgeLess(Parc a,Parc b)
160 std::vector<Edge> arcs;
162 namespace _delaunay_bits {
165 int prev, curr, next;
167 Part(int p, int c, int n) : prev(p), curr(c), next(n) {}
170 inline std::ostream& operator<<(std::ostream& os, const Part& part) {
171 os << '(' << part.prev << ',' << part.curr << ',' << part.next << ')';
175 inline double circle_point(const Point& p, const Point& q, const Point& r) {
176 double a = p.x * (q.y - r.y) + q.x * (r.y - p.y) + r.x * (p.y - q.y);
177 if (a == 0) return std::numeric_limits<double>::quiet_NaN();
179 double d = (p.x * p.x + p.y * p.y) * (q.y - r.y) +
180 (q.x * q.x + q.y * q.y) * (r.y - p.y) +
181 (r.x * r.x + r.y * r.y) * (p.y - q.y);
183 double e = (p.x * p.x + p.y * p.y) * (q.x - r.x) +
184 (q.x * q.x + q.y * q.y) * (r.x - p.x) +
185 (r.x * r.x + r.y * r.y) * (p.x - q.x);
187 double f = (p.x * p.x + p.y * p.y) * (q.x * r.y - r.x * q.y) +
188 (q.x * q.x + q.y * q.y) * (r.x * p.y - p.x * r.y) +
189 (r.x * r.x + r.y * r.y) * (p.x * q.y - q.x * p.y);
191 return d / (2 * a) + sqrt((d * d + e * e) / (4 * a * a) + f / a);
194 inline bool circle_form(const Point& p, const Point& q, const Point& r) {
195 return rot90(q - p) * (r - q) < 0.0;
198 inline double intersection(const Point& p, const Point& q, double sx) {
199 const double epsilon = 1e-8;
201 if (p.x == q.x) return (p.y + q.y) / 2.0;
203 if (sx < p.x + epsilon) return p.y;
204 if (sx < q.x + epsilon) return q.y;
206 double a = q.x - p.x;
207 double b = (q.x - sx) * p.y - (p.x - sx) * q.y;
208 double d = (q.x - sx) * (p.x - sx) * (p - q).normSquare();
209 return (b - sqrt(d)) / a;
215 YLess(const std::vector<Point>& points, double& sweep)
216 : _points(points), _sweep(sweep) {}
218 bool operator()(const Part& l, const Part& r) const {
219 const double epsilon = 1e-8;
221 // std::cerr << l << " vs " << r << std::endl;
222 double lbx = l.prev != -1 ?
223 intersection(_points[l.prev], _points[l.curr], _sweep) :
224 - std::numeric_limits<double>::infinity();
225 double rbx = r.prev != -1 ?
226 intersection(_points[r.prev], _points[r.curr], _sweep) :
227 - std::numeric_limits<double>::infinity();
228 double lex = l.next != -1 ?
229 intersection(_points[l.curr], _points[l.next], _sweep) :
230 std::numeric_limits<double>::infinity();
231 double rex = r.next != -1 ?
232 intersection(_points[r.curr], _points[r.next], _sweep) :
233 std::numeric_limits<double>::infinity();
235 if (lbx > lex) std::swap(lbx, lex);
236 if (rbx > rex) std::swap(rbx, rex);
238 if (lex < epsilon + rex && lbx + epsilon < rex) return true;
239 if (rex < epsilon + lex && rbx + epsilon < lex) return false;
243 const std::vector<Point>& _points;
249 typedef std::multimap<double, BeachIt> SpikeHeap;
251 typedef std::multimap<Part, SpikeHeap::iterator, YLess> Beach;
256 BeachIt(Beach::iterator iter) : it(iter) {}
261 inline void delaunay() {
262 Counter cnt("Number of arcs added: ");
264 using namespace _delaunay_bits;
266 typedef _delaunay_bits::Part Part;
267 typedef std::vector<std::pair<double, int> > SiteHeap;
270 std::vector<Point> points;
271 std::vector<Node> nodes;
273 for (NodeIt it(g); it != INVALID; ++it) {
275 points.push_back(coords[it]);
278 SiteHeap siteheap(points.size());
283 for (int i = 0; i < int(siteheap.size()); ++i) {
284 siteheap[i] = std::make_pair(points[i].x, i);
287 std::sort(siteheap.begin(), siteheap.end());
288 sweep = siteheap.front().first;
290 YLess yless(points, sweep);
295 std::set<std::pair<int, int> > arcs;
301 while (siteindex < int(siteheap.size()) &&
302 siteheap[0].first == siteheap[siteindex].first) {
303 front.push_back(std::make_pair(points[siteheap[siteindex].second].y,
304 siteheap[siteindex].second));
308 std::sort(front.begin(), front.end());
310 for (int i = 0; i < int(front.size()); ++i) {
311 int prev = (i == 0 ? -1 : front[i - 1].second);
312 int curr = front[i].second;
313 int next = (i + 1 == int(front.size()) ? -1 : front[i + 1].second);
315 beach.insert(std::make_pair(Part(prev, curr, next),
320 while (siteindex < int(points.size()) || !spikeheap.empty()) {
322 SpikeHeap::iterator spit = spikeheap.begin();
324 if (siteindex < int(points.size()) &&
325 (spit == spikeheap.end() || siteheap[siteindex].first < spit->first)) {
326 int site = siteheap[siteindex].second;
327 sweep = siteheap[siteindex].first;
329 Beach::iterator bit = beach.upper_bound(Part(site, site, site));
331 if (bit->second != spikeheap.end()) {
332 spikeheap.erase(bit->second);
335 int prev = bit->first.prev;
336 int curr = bit->first.curr;
337 int next = bit->first.next;
341 SpikeHeap::iterator pit = spikeheap.end();
343 circle_form(points[prev], points[curr], points[site])) {
344 double x = circle_point(points[prev], points[curr], points[site]);
345 pit = spikeheap.insert(std::make_pair(x, BeachIt(beach.end())));
347 beach.insert(std::make_pair(Part(prev, curr, site), pit));
349 beach.insert(std::make_pair(Part(prev, curr, site), pit));
352 beach.insert(std::make_pair(Part(curr, site, curr), spikeheap.end()));
354 SpikeHeap::iterator nit = spikeheap.end();
356 circle_form(points[site], points[curr],points[next])) {
357 double x = circle_point(points[site], points[curr], points[next]);
358 nit = spikeheap.insert(std::make_pair(x, BeachIt(beach.end())));
360 beach.insert(std::make_pair(Part(site, curr, next), nit));
362 beach.insert(std::make_pair(Part(site, curr, next), nit));
369 Beach::iterator bit = spit->second.it;
371 int prev = bit->first.prev;
372 int curr = bit->first.curr;
373 int next = bit->first.next;
376 std::pair<int, int> arc;
379 std::make_pair(prev, curr) : std::make_pair(curr, prev);
381 if (arcs.find(arc) == arcs.end()) {
383 g.addEdge(nodes[prev], nodes[curr]);
388 std::make_pair(curr, next) : std::make_pair(next, curr);
390 if (arcs.find(arc) == arcs.end()) {
392 g.addEdge(nodes[curr], nodes[next]);
397 Beach::iterator pbit = bit; --pbit;
398 int ppv = pbit->first.prev;
399 Beach::iterator nbit = bit; ++nbit;
400 int nnt = nbit->first.next;
402 if (bit->second != spikeheap.end()) spikeheap.erase(bit->second);
403 if (pbit->second != spikeheap.end()) spikeheap.erase(pbit->second);
404 if (nbit->second != spikeheap.end()) spikeheap.erase(nbit->second);
410 SpikeHeap::iterator pit = spikeheap.end();
411 if (ppv != -1 && ppv != next &&
412 circle_form(points[ppv], points[prev], points[next])) {
413 double x = circle_point(points[ppv], points[prev], points[next]);
414 if (x < sweep) x = sweep;
415 pit = spikeheap.insert(std::make_pair(x, BeachIt(beach.end())));
417 beach.insert(std::make_pair(Part(ppv, prev, next), pit));
419 beach.insert(std::make_pair(Part(ppv, prev, next), pit));
422 SpikeHeap::iterator nit = spikeheap.end();
423 if (nnt != -1 && prev != nnt &&
424 circle_form(points[prev], points[next], points[nnt])) {
425 double x = circle_point(points[prev], points[next], points[nnt]);
426 if (x < sweep) x = sweep;
427 nit = spikeheap.insert(std::make_pair(x, BeachIt(beach.end())));
429 beach.insert(std::make_pair(Part(prev, next, nnt), nit));
431 beach.insert(std::make_pair(Part(prev, next, nnt), nit));
437 for (Beach::iterator it = beach.begin(); it != beach.end(); ++it) {
438 int curr = it->first.curr;
439 int next = it->first.next;
441 if (next == -1) continue;
443 std::pair<int, int> arc;
446 std::make_pair(curr, next) : std::make_pair(next, curr);
448 if (arcs.find(arc) == arcs.end()) {
450 g.addEdge(nodes[curr], nodes[next]);
458 Counter cnt("Number of arcs removed: ");
459 Bfs<ListGraph> bfs(g);
460 for(std::vector<Edge>::reverse_iterator ei=arcs.rbegin();
461 ei!=arcs.rend();++ei)
467 if(bfs.predArc(b)==INVALID || bfs.dist(b)>d)
475 Counter cnt("Number of arcs removed: ");
476 for(std::vector<Edge>::reverse_iterator ei=arcs.rbegin();
477 ei!=arcs.rend();++ei)
482 ConstMap<Arc,int> cegy(1);
483 Suurballe<ListGraph,ConstMap<Arc,int> > sur(g,cegy,a,b);
485 if(k<2 || sur.totalLength()>d)
488 // else std::cout << "Remove arc " << g.id(a) << "-" << g.id(b) << '\n';
492 void sparseTriangle(int d)
494 Counter cnt("Number of arcs added: ");
495 std::vector<Parc> pedges;
496 for(NodeIt n(g);n!=INVALID;++n)
497 for(NodeIt m=++(NodeIt(n));m!=INVALID;++m)
502 p.len=(coords[m]-coords[n]).normSquare();
505 std::sort(pedges.begin(),pedges.end(),pedgeLess);
506 for(std::vector<Parc>::iterator pi=pedges.begin();pi!=pedges.end();++pi)
508 Line li(pi->a,pi->b);
510 for(;e!=INVALID && !cross(e,li);++e) ;
513 ConstMap<Arc,int> cegy(1);
514 Suurballe<ListGraph,ConstMap<Arc,int> >
515 sur(g,cegy,pi->a,pi->b);
517 if(k<2 || sur.totalLength()>d)
519 ne=g.addEdge(pi->a,pi->b);
527 template <typename Graph, typename CoordMap>
528 class LengthSquareMap {
530 typedef typename Graph::Edge Key;
531 typedef typename CoordMap::Value::Value Value;
533 LengthSquareMap(const Graph& graph, const CoordMap& coords)
534 : _graph(graph), _coords(coords) {}
536 Value operator[](const Key& key) const {
537 return (_coords[_graph.v(key)] -
538 _coords[_graph.u(key)]).normSquare();
544 const CoordMap& _coords;
548 std::vector<Parc> pedges;
550 std::cout << T.realTime() << "s: Creating delaunay triangulation...\n";
552 std::cout << T.realTime() << "s: Calculating spanning tree...\n";
553 LengthSquareMap<ListGraph, ListGraph::NodeMap<Point> > ls(g, coords);
554 ListGraph::EdgeMap<bool> tree(g);
555 kruskal(g, ls, tree);
556 std::cout << T.realTime() << "s: Removing non tree arcs...\n";
557 std::vector<Edge> remove;
558 for (EdgeIt e(g); e != INVALID; ++e) {
559 if (!tree[e]) remove.push_back(e);
561 for(int i = 0; i < int(remove.size()); ++i) {
564 std::cout << T.realTime() << "s: Done\n";
569 std::cout << "Find a tree..." << std::endl;
573 std::cout << "Total arc length (tree) : " << totalLen() << std::endl;
575 std::cout << "Make it Euler..." << std::endl;
578 std::vector<Node> leafs;
579 for(NodeIt n(g);n!=INVALID;++n)
580 if(countIncEdges(g,n)%2==1) leafs.push_back(n);
582 // for(unsigned int i=0;i<leafs.size();i+=2)
583 // g.addArc(leafs[i],leafs[i+1]);
585 std::vector<Parc> pedges;
586 for(unsigned int i=0;i<leafs.size()-1;i++)
587 for(unsigned int j=i+1;j<leafs.size();j++)
594 p.len=(coords[m]-coords[n]).normSquare();
597 std::sort(pedges.begin(),pedges.end(),pedgeLess);
598 for(unsigned int i=0;i<pedges.size();i++)
599 if(countIncEdges(g,pedges[i].a)%2 &&
600 countIncEdges(g,pedges[i].b)%2)
601 g.addEdge(pedges[i].a,pedges[i].b);
604 for(NodeIt n(g);n!=INVALID;++n)
605 if(countIncEdges(g,n)%2 || countIncEdges(g,n)==0 )
606 std::cout << "GEBASZ!!!" << std::endl;
608 for(EdgeIt e(g);e!=INVALID;++e)
610 std::cout << "LOOP GEBASZ!!!" << std::endl;
612 std::cout << "Number of arcs : " << countEdges(g) << std::endl;
614 std::cout << "Total arc length (euler) : " << totalLen() << std::endl;
616 ListGraph::EdgeMap<Arc> enext(g);
618 EulerIt<ListGraph> e(g);
621 // std::cout << "Tour arc: " << g.id(Edge(e)) << std::endl;
622 for(++e;e!=INVALID;++e)
624 // std::cout << "Tour arc: " << g.id(Edge(e)) << std::endl;
631 std::cout << "Creating a tour from that..." << std::endl;
633 int nnum = countNodes(g);
634 int ednum = countEdges(g);
636 for(Arc p=enext[EdgeIt(g)];ednum>nnum;p=enext[p])
638 // std::cout << "Checking arc " << g.id(p) << std::endl;
642 Node n1=g.oppositeNode(n2,e);
643 Node n3=g.oppositeNode(n2,f);
644 if(countIncEdges(g,n2)>2)
646 // std::cout << "Remove an Arc" << std::endl;
652 Arc ne=g.direct(g.addEdge(n1,n3),n1);
664 std::cout << "Total arc length (tour) : " << totalLen() << std::endl;
666 std::cout << "2-opt the tour..." << std::endl;
670 std::cout << "Total arc length (2-opt tour) : " << totalLen() << std::endl;
674 int main(int argc,const char **argv)
676 ArgParser ap(argc,argv);
679 bool disc_d, square_d, gauss_d;
680 // bool tsp_a,two_a,tree_a;
685 std::string ndist("disc");
686 ap.refOption("n", "Number of nodes (default is 100)", N)
687 .intOption("g", "Girth parameter (default is 10)", 10)
688 .refOption("cities", "Number of cities (default is 1)", num_of_cities)
689 .refOption("area", "Full relative area of the cities (default is 1)", area)
690 .refOption("disc", "Nodes are evenly distributed on a unit disc (default)",disc_d)
691 .optionGroup("dist", "disc")
692 .refOption("square", "Nodes are evenly distributed on a unit square", square_d)
693 .optionGroup("dist", "square")
695 "Nodes are located according to a two-dim gauss distribution",
697 .optionGroup("dist", "gauss")
698 // .mandatoryGroup("dist")
699 .onlyOneGroup("dist")
700 .boolOption("eps", "Also generate .eps output (prefix.eps)")
701 .boolOption("nonodes", "Draw the edges only in the generated .eps")
702 .boolOption("dir", "Directed digraph is generated (each arcs are replaced by two directed ones)")
703 .boolOption("2con", "Create a two connected planar digraph")
704 .optionGroup("alg","2con")
705 .boolOption("tree", "Create a min. cost spanning tree")
706 .optionGroup("alg","tree")
707 .boolOption("tsp", "Create a TSP tour")
708 .optionGroup("alg","tsp")
709 .boolOption("tsp2", "Create a TSP tour (tree based)")
710 .optionGroup("alg","tsp2")
711 .boolOption("dela", "Delaunay triangulation digraph")
712 .optionGroup("alg","dela")
714 .boolOption("rand", "Use time seed for random number generator")
715 .optionGroup("rand", "rand")
716 .intOption("seed", "Random seed", -1)
717 .optionGroup("rand", "seed")
718 .onlyOneGroup("rand")
719 .other("[prefix]","Prefix of the output files. Default is 'lgf-gen-out'")
724 std::cout << "Random number seed: " << seed << std::endl;
727 if (ap.given("seed")) {
728 int seed = ap["seed"];
729 std::cout << "Random number seed: " << seed << std::endl;
734 switch(ap.files().size())
737 prefix="lgf-gen-out";
740 prefix=ap.files()[0];
743 std::cerr << "\nAt most one prefix can be given\n\n";
748 std::vector<double> sizes;
749 std::vector<double> cum_sizes;
750 for(int s=0;s<num_of_cities;s++)
752 // sum_sizes+=rnd.exponential();
756 cum_sizes.push_back(sum_sizes);
759 for(int s=0;s<num_of_cities;s++)
761 Point center=(num_of_cities==1?Point(0,0):rnd.disc());
763 for(;i<N*(cum_sizes[s]/sum_sizes);i++) {
766 coords[n]=center+rnd.gauss2()*area*
767 std::sqrt(sizes[s]/sum_sizes);
770 for(;i<N*(cum_sizes[s]/sum_sizes);i++) {
773 coords[n]=center+Point(rnd()*2-1,rnd()*2-1)*area*
774 std::sqrt(sizes[s]/sum_sizes);
776 else if(disc_d || true)
777 for(;i<N*(cum_sizes[s]/sum_sizes);i++) {
780 coords[n]=center+rnd.disc()*area*
781 std::sqrt(sizes[s]/sum_sizes);
785 // for (ListGraph::NodeIt n(g); n != INVALID; ++n) {
786 // std::cerr << coords[n] << std::endl;
791 std::cout << "#2-opt improvements: " << tsp_impr_num << std::endl;
795 std::cout << "#2-opt improvements: " << tsp_impr_num << std::endl;
797 else if(ap["2con"]) {
798 std::cout << "Make triangles\n";
800 sparseTriangle(ap["g"]);
801 std::cout << "Make it sparser\n";
804 else if(ap["tree"]) {
807 else if(ap["dela"]) {
812 std::cout << "Number of nodes : " << countNodes(g) << std::endl;
813 std::cout << "Number of arcs : " << countEdges(g) << std::endl;
815 for(EdgeIt e(g);e!=INVALID;++e)
816 tlen+=sqrt((coords[g.v(e)]-coords[g.u(e)]).normSquare());
817 std::cout << "Total arc length : " << tlen << std::endl;
820 graphToEps(g,prefix+".eps").scaleToA4().
821 scale(600).nodeScale(.005).arcWidthScale(.001).preScale(false).
822 coords(coords).hideNodes(ap.given("nonodes")).run();
825 DigraphWriter<ListGraph>(g,prefix+".lgf").
826 nodeMap("coordinates_x",scaleMap(xMap(coords),600)).
827 nodeMap("coordinates_y",scaleMap(yMap(coords),600)).
829 else GraphWriter<ListGraph>(g,prefix+".lgf").
830 nodeMap("coordinates_x",scaleMap(xMap(coords),600)).
831 nodeMap("coordinates_y",scaleMap(yMap(coords),600)).