File: mst.cpp

package info (click to toggle)
cgal 6.1.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 144,952 kB
  • sloc: cpp: 811,597; ansic: 208,576; sh: 493; python: 411; makefile: 286; javascript: 174
file content (109 lines) | stat: -rw-r--r-- 3,554 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
// Copyright (c) 2018 Inria
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org)
//
// SPDX-License-Identifier: LGPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s)     : Marc Glisse

#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Delaunay_triangulation_2.h>
#include <CGAL/CGAL_Ipelet_base.h>
#include <CGAL/boost/graph/graph_traits_Delaunay_triangulation_2.h>
#include <boost/graph/kruskal_min_spanning_tree.hpp>
#include <boost/graph/filtered_graph.hpp>

namespace CGAL_mst {

typedef CGAL::Exact_predicates_inexact_constructions_kernel Kernel;
typedef CGAL::Delaunay_triangulation_2<Kernel>              Triangulation;

const std::string Slab[] = {
  "MST", "Help"
};

const std::string Hmsg[] = {
  "Draw the minimum spanning tree of a set of points"
};

struct mstIpelet
  : CGAL::Ipelet_base<Kernel, 2> {
  mstIpelet() : CGAL::Ipelet_base<Kernel, 2>("Minimum spanning tree", Slab, Hmsg){}
  void protected_run(int);
};

template <typename T>
struct Is_finite {
  const T* t_;
  Is_finite()
    : t_(nullptr)
  {}
  Is_finite(const T& t)
    : t_(&t)
  { }
  template <typename VertexOrEdge>
  bool operator()(const VertexOrEdge& voe) const {
    return ! t_->is_infinite(voe);
  }
};
typedef Is_finite<Triangulation> Filter;
typedef boost::filtered_graph<Triangulation,Filter,Filter> Finite_triangulation;
typedef boost::graph_traits<Finite_triangulation>::vertex_descriptor vertex_descriptor;
typedef boost::graph_traits<Finite_triangulation>::vertex_iterator vertex_iterator;
typedef boost::graph_traits<Finite_triangulation>::edge_descriptor edge_descriptor;
// The BGL makes use of indices associated to the vertices
// We use a std::map to store the index
typedef std::map<vertex_descriptor,int> VertexIndexMap;
VertexIndexMap vertex_id_map;
// A std::map is not a property map, because it is not lightweight
typedef boost::associative_property_map<VertexIndexMap> VertexIdPropertyMap;
VertexIdPropertyMap vertex_index_pmap(vertex_id_map);

void mstIpelet::protected_run(int /*fn*/)
{
  std::list<Point_2> pt_list;

  read_active_objects(
    CGAL::dispatch_or_drop_output<Point_2>(
      std::back_inserter(pt_list)
    )
  );

  if (pt_list.empty()) {
    print_error_message("No mark selected");
    return;
  }

  Triangulation t(pt_list.begin(), pt_list.end());
  Filter is_finite(t);
  Finite_triangulation ft(t, is_finite, is_finite);

  vertex_iterator vit, ve;
  // Associate indices to the vertices
  int index = 0;
  // std::tie assigns the first and second element of the std::pair
  // returned by boost::vertices to the variables vit and ve
  for(std::tie(vit,ve)=boost::vertices(ft); vit!=ve; ++vit ){
    vertex_descriptor vd = *vit;
    vertex_id_map[vd] = index++;
    }
  // We use the default edge weight which is the squared length of the edge
  // This property map is defined in graph_traits_Triangulation_2.h
  // In the function call you can see a named parameter: vertex_index_map
   std::list<edge_descriptor> mst;
   boost::kruskal_minimum_spanning_tree(ft,
                    std::back_inserter(mst),
                    vertex_index_map(vertex_index_pmap));
   for(std::list<edge_descriptor>::iterator it = mst.begin(); it != mst.end(); ++it){
     edge_descriptor ed = *it;
     vertex_descriptor svd = source(ed,t);
     vertex_descriptor tvd = target(ed,t);
     Triangulation::Vertex_handle sv = svd;
     Triangulation::Vertex_handle tv = tvd;
     draw_in_ipe(Kernel::Segment_2(sv->point(), tv->point()));
   }
}
}

CGAL_IPELET(CGAL_mst::mstIpelet)