File: corefinement_mesh_union_progress.cpp

package info (click to toggle)
cgal 6.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 144,912 kB
  • sloc: cpp: 810,858; ansic: 208,477; sh: 493; python: 411; makefile: 286; javascript: 174
file content (219 lines) | stat: -rw-r--r-- 5,646 bytes parent folder | download | duplicates (2)
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
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Surface_mesh.h>

#include <CGAL/Polygon_mesh_processing/corefinement.h>
#include <CGAL/Polygon_mesh_processing/IO/polygon_mesh_io.h>

#include <fstream>

typedef CGAL::Exact_predicates_inexact_constructions_kernel   K;
typedef CGAL::Surface_mesh<K::Point_3>                        Mesh;

namespace PMP = CGAL::Polygon_mesh_processing;

struct Visitor_rep{

  Visitor_rep(double normalize = 4)
    : normalize(normalize)
  {
    t.start();
  }

  void progress_filtering_intersections(double d)
   {
    d /= normalize;
    total += d;
    if(total > bound){
      std::cout << std::setprecision(3) << total*100 << " %   in " << std::setprecision(5) << t.time() << " sec." << std::endl;
      bound += 0.1;
    }
  }

  void start_triangulating_faces(std::size_t tf)
  {
    tfaces = tf;
    bound_faces = tf/10;
  }

  void face_triangulation(std::size_t i)
  {
    if(i> bound_faces){
      std::cout << double(i)/double(tfaces) * 100  << " %" << std::endl;
      bound_faces += tfaces/10;
    }
  }

  void start_coplanar_faces(std::size_t tc)
  {
    std::cout << "Visitor::start_coplanar_faces() at " << t.time() << " sec." << std::endl;
    tcoplanar= tc;
    count_coplanar = 0;
    bound_coplanar = tcoplanar/10;
  }

  void intersection_of_coplanar_faces_step()
  {
    ++count_coplanar;
    if(count_coplanar> bound_coplanar){
      std::cout << "Visitor::coplanar_faces: " << double(count_coplanar)/double(tcoplanar) * 100  << " % " << std::endl;
      bound_coplanar += tcoplanar/10;
    }
  }

  void start_intersection_points(std::size_t ti)
  {
    std::cout << "Visitor::start_intersection_points() at " << t.time() << " sec." << std::endl;
    tintersection= ti;
    count_intersection = 0;
    bound_intersection = tintersection/10;
  }

  void edge_face_intersections_step()
  {
    ++count_intersection;
    if(count_intersection> bound_intersection){
      std::cout << "Visitor::intersection_points: " << double(count_intersection)/double(tintersection) * 100  << " % " << std::endl;
      bound_intersection += tintersection/10;
    }
  }

  double time() const
  {
    return t.time();
  }

  double normalize;
  double bound = 0.1;
  double total = 0;
  std::size_t count = 0;

  std::size_t bound_faces = 0;
  std::size_t tfaces = 0;

  std::size_t bound_coplanar = 0;
  std::size_t tcoplanar = 0;
  std::size_t count_coplanar = 0;

  std::size_t bound_intersection = 0;
  std::size_t tintersection = 0;
  std::size_t count_intersection = 0;
  CGAL::Timer t;
};


struct Visitor :
  public PMP::Corefinement::Default_visitor<Mesh>
{
  std::shared_ptr<Visitor_rep> sptr;
  mutable std::size_t tf_counter = 0;

  Visitor()
    : sptr(std::make_shared<Visitor_rep>())
  {}

  void progress_filtering_intersections(double d)
  {
    sptr->progress_filtering_intersections(d);
  }

  void start_filtering_intersections() const
  {
    std::cout << "Visitor::start_filtering_intersections() at " << sptr->time() << " sec." << std::endl;
  }
  void end_filtering_intersections() const
  {
    std::cout << "Visitor::end_filtering_intersections() at " << sptr->time() << " sec."  << std::endl;
  }

  void start_triangulating_faces(std::size_t tf) const
  {
    std::cout << "Visitor::start_triangulation() with " << tf << " faces at " << sptr->time() << " sec."  << std::endl;
    sptr->start_triangulating_faces(tf);
    tf_counter = 0;
  }

  void triangulating_faces_step() const
  {
    sptr->face_triangulation(tf_counter++);
  }

  void end_triangulating_faces()const
  {
    std::cout << "Visitor::end_triangulating_faces() at " << sptr->time() << " sec."  << std::endl;
  }

  void start_handling_intersection_of_coplanar_faces(std::size_t i) const
  {
    sptr->start_coplanar_faces(i);
  }

  void intersection_of_coplanar_faces_step() const
  {
    sptr->intersection_of_coplanar_faces_step();
  }

  void end_handling_intersection_of_coplanar_faces() const
  {
    std::cout << "Visitor::end_coplanar_faces() at " << sptr->time() << " sec." << std::endl;
  }

  void start_handling_edge_face_intersections(std::size_t i) const
  {
    sptr->start_intersection_points(i);
  }

  void edge_face_intersections_step() const
  {
    sptr->edge_face_intersections_step();
  }

  void end_handling_edge_face_intersections() const
  {
    std::cout << "Visitor::end_intersection_points() at " << sptr->time() << " sec." << std::endl;
  }

  void start_building_output() const
  {
    std::cout << "Visitor::start_building_output() at " << sptr->time()  << " sec."<< std::endl;
  }

  void end_building_output() const
  {
    std::cout << "Visitor::end_building_output() at " << sptr->time() << " sec." << std::endl;
  }
};


int main(int argc, char* argv[])
{
  const std::string filename1 = (argc > 1) ? argv[1] : CGAL::data_file_path("meshes/blobby.off");
  const std::string filename2 = (argc > 2) ? argv[2] : CGAL::data_file_path("meshes/eight.off");

  Mesh mesh1, mesh2;
  if(!CGAL::IO::read_polygon_mesh(filename1, mesh1) || !CGAL::IO::read_polygon_mesh(filename2, mesh2))
  {
    std::cerr << "Invalid input." << std::endl;
    return 1;
  }

  CGAL::Timer t;
  t.start();
  Mesh out;
  Visitor visitor;

  bool valid_union = PMP::corefine_and_compute_union (mesh1,mesh2, out, CGAL::parameters::visitor(visitor));

  std::cout << "Global timer = " << t.time() << " sec." << std::endl;


  if(valid_union)
  {
    std::cout << "Union was successfully computed\n";
    CGAL::IO::write_polygon_mesh("union.off", out, CGAL::parameters::stream_precision(17));
    return 0;
  }

  std::cout << "Union could not be computed\n";

  return 1;
}