MEPP2 Project
cmdm.h
Go to the documentation of this file.
1 // Copyright (c) 2012-2020 University of Lyon and CNRS (France).
2 // All rights reserved.
3 //
4 // This file is part of MEPP2; you can redistribute it and/or modify
5 // it under the terms of the GNU General Public License as published
6 // by the Free Software Foundation; either version 3 of the License,
7 // or (at your option) any later version.
8 //
9 // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE
10 // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE.
11 
12 #pragma once
13 
17 
18 #include <CGAL/AABB_tree.h>
19 #include <CGAL/AABB_traits.h>
20 #include <CGAL/AABB_face_graph_triangle_primitive.h>
21 
24 #include "curvature_cmdm.hpp"
28 
29 #include <boost/foreach.hpp>
30 
31 #ifdef _MSC_VER
32 // disable some warnings on Windows
33 #pragma warning(push)
34 #pragma warning(disable : 4244)
35 #pragma warning(disable : 4267)
36 // 4244 & 4267: converting type A to type B, possible data loss
37 #endif
38 
39 namespace FEVV {
40 
41 using Triangle = CGALKernel::Triangle_3;
42 using Segment = CGALKernel::Segment_3;
43 using Ray = CGALKernel::Ray_3;
44 using Vertex_Descriptor = MeshSurface::Vertex_index;
45 
46 namespace Filters {
47 
48 template< typename Point3d, typename Vector >
49 bool
50 sphere_clip_vector_cmdm(Point3d &o, double r, const Point3d &p, Vector &v)
51 {
52 
53  Vector w = p - o;
54  double a = (v * v);
55  double b = 2.0 * v * w;
56  double c = (w * w) - r * r;
57  double delta = b * b - 4 * a * c;
58 
59  if(a == 0)
60  return true;
61 
62  if(delta < 0)
63  {
64  // Should not happen, but happens sometimes (numerical precision)
65 
66  return true;
67  }
68  double t = (-b + std::sqrt(delta)) / (2.0 * a);
69  if(t < 0.0)
70  {
71 
72  // Should not happen, but happens sometimes (numerical precision)
73  return true;
74  }
75  if(t >= 1.0)
76  {
77  // Inside the sphere
78  return false;
79  }
80 
81  if(t < 1.e-16)
82  {
83 
84  t = 0.01;
85  }
86 
87  v = v * t;
88 
89  return true;
90 }
91 
107 template< typename VertexIterator,
108  typename HalfedgeGraph,
109  typename PointMap,
110  typename CMDMNearestMap >
111 void
112 compute_geometry_statistics(const HalfedgeGraph &/*m_degr*/,
113  const PointMap &pm_degr,
114  const CMDMNearestMap &nearest_pmap,
115  const VertexIterator p_vertex,
116  std::vector< double > &tab_distance1,
117  std::vector< double > &tab_distance2,
118  std::vector< CGALPoint > &tab_point1,
119  std::vector< CGALPoint > &tab_point2,
120  double *tab_f1,
121  double *tab_f2,
122  double *tab_f3,
123  double radius)
124 {
125  double moyenne1 = 0;
126  double variance1 = 0;
127 
128  double moyenne2 = 0;
129  double variance2 = 0;
130 
131  double covariance = 0;
132  double variance = radius / 3;
133  std::vector< double > tab_wi1;
134  std::vector< double > tab_wi2;
135  double somme_distance1 = 0;
136  double somme_distance2 = 0;
137  double somme_wi1 = 0;
138  double somme_wi2 = 0;
139  for(unsigned int i = 0; i < tab_point1.size(); i++)
140  {
141  auto distance_pt1 = tab_point1[i] - get(pm_degr, p_vertex);
142  double dist_pt1 = sqrt(distance_pt1 * distance_pt1);
143  double wi1 = 1 / variance / sqrt(2 * 3.141592) *
144  exp(-(dist_pt1 * dist_pt1) / 2 / variance / variance);
145 
146  tab_wi1.push_back(wi1);
147  somme_wi1 += wi1;
148  somme_distance1 += tab_distance1[i] * wi1;
149 
150  auto distance_pt2 = tab_point2[i] - get(nearest_pmap, p_vertex);
151  double dist_pt2 = sqrt(distance_pt2 * distance_pt2);
152 
153  double wi2 = 1 / variance / sqrt(2 * 3.141592) *
154  exp(-(dist_pt2 * dist_pt2) / 2 / variance / variance);
155  tab_wi2.push_back(wi2);
156  somme_wi2 += wi2;
157  somme_distance2 += tab_distance2[i] * wi2;
158  }
159  moyenne1 = somme_distance1 / (double)somme_wi1;
160  moyenne2 = somme_distance2 / (double)somme_wi2;
161 
162 
163  for(unsigned int i = 0; i < tab_point1.size(); i++)
164  {
165  variance1 += tab_wi1[i] * pow(tab_distance1[i] - moyenne1, 2);
166  variance2 += tab_wi2[i] * pow(tab_distance2[i] - moyenne2, 2);
167 
168  covariance += tab_wi1[i] * (tab_distance1[i] - moyenne1) *
169  (tab_distance2[i] - moyenne2);
170  }
171 
172  variance1 = variance1 / somme_wi1;
173  variance2 = variance2 / somme_wi2;
174  variance1 = sqrt(variance1);
175  variance2 = sqrt(variance2);
176  covariance = covariance / somme_wi1;
177  covariance = (covariance < 0) ? 0 : covariance;
178 
179  // We then compute the local geometry-based features
180  // 1) Curvature comparison
181  double f1 =
182  (std::fabs(moyenne1 - moyenne2)) / (std::max(moyenne1, moyenne2) + 1);
183  // 2) Curvature contrast
184  double f2 =
185  (std::fabs(variance1 - variance2)) / (std::max(variance1, variance2) + 1);
186  // 3) Curvature structure
187  double f3 = (std::fabs(variance1 * variance2 - covariance)) /
188  (variance1 * variance2 + 1);
189 
190  tab_f1[p_vertex] += f1;
191  tab_f2[p_vertex] += f2;
192  tab_f3[p_vertex] += f3;
193 }
194 
209 template< typename VertexIterator,
210  typename HalfedgeGraph,
211  typename PointMap,
212  typename CMDMNearestMap >
213 void
214 compute_color_statistics(const HalfedgeGraph &/*m_degr*/,
215  const PointMap &pm_degr,
216  const CMDMNearestMap &nearest_pmap,
217  const VertexIterator p_vertex,
218  std::vector< std::array< double, 3 > > &tab_color1,
219  std::vector< std::array< double, 3 > > &tab_color2,
220  std::vector< CGALPoint > &tab_point1,
221  std::vector< CGALPoint > &tab_point2,
222  double *tab_f4,
223  double *tab_f5,
224  double *tab_f6,
225  double *tab_f7,
226  double *tab_f8,
227  double radius)
228 {
229  double variance = radius / 3;
230  std::vector< double > tab_wi1;
231  std::vector< double > tab_wi2;
232  double somme_wi1 = 0;
233  double somme_wi2 = 0;
234 
235  for(unsigned int i = 0; i < tab_point1.size(); i++)
236  {
237  auto distance_pt1 = tab_point1[i] - get(pm_degr, p_vertex);
238  double dist_pt1 = sqrt(distance_pt1 * distance_pt1);
239  double wi1 = 1 / variance / sqrt(2 * 3.141592) *
240  exp(-(dist_pt1 * dist_pt1) / 2 / variance / variance);
241  tab_wi1.push_back(wi1);
242  somme_wi1 += wi1;
243 
244  auto distance_pt2 = tab_point2[i] - get(nearest_pmap, p_vertex);
245  double dist_pt2 = sqrt(distance_pt2 * distance_pt2);
246  double wi2 = 1 / variance / sqrt(2 * 3.141592) *
247  exp(-(dist_pt2 * dist_pt2) / 2 / variance / variance);
248  tab_wi2.push_back(wi2);
249  somme_wi2 += wi2;
250  }
251 
252  std::vector< double > tab_Chr1; // Chroma from the neighborhood of pVertex
253  // regarding the first polyhedron
254  std::vector< double > tab_Chr2; // Chroma from the neighborhood of pVertex
255  // regarding the second polyhedron
256  std::vector< double > tab_dH; // Hue difference
257 
258  double muL1 = 0, muCh1 = 0;
259  double variance1_L = 0, sL1 = 0;
260 
261  double muL2 = 0, muCh2 = 0;
262  double variance2_L = 0, sL2 = 0;
263 
264  double dL_sq = 0, dCh_sq = 0, dH_sq = 0, sL12 = 0; // Mixed terms
265 
266  double somme1_L = 0, somme1_Chr = 0;
267  double somme2_L = 0, somme2_Chr = 0;
268  for(unsigned int i = 0; i < tab_point1.size(); i++)
269  {
270  // Chroma = sqrt(A^2 + B^2)
271  tab_Chr1.push_back(sqrt(tab_color1[i][1] * tab_color1[i][1] +
272  tab_color1[i][2] * tab_color1[i][2]));
273  somme1_L += tab_color1[i][0] * tab_wi1[i];
274  somme1_Chr += tab_Chr1[i] * tab_wi1[i];
275 
276  tab_Chr2.push_back(sqrt(tab_color2[i][1] * tab_color2[i][1] +
277  tab_color2[i][2] * tab_color2[i][2]));
278  somme2_L += tab_color2[i][0] * tab_wi2[i];
279  somme2_Chr += tab_Chr2[i] * tab_wi2[i];
280 
281  // Hue difference = sqrt[(A1-A2)^2 + (B1-B2)^2 - (Chr1-Chr2)^2]
282  double dH = 0;
283  dH = std::pow(tab_color1[i][1] - tab_color2[i][1], 2) +
284  std::pow(tab_color1[i][2] - tab_color2[i][2], 2) -
285  std::pow(tab_Chr1[i] - tab_Chr2[i], 2);
286  tab_dH.push_back((dH < 0) ? 0 : sqrt(dH));
287  }
288 
289  muL1 = somme1_L / somme_wi1;
290  muCh1 = somme1_Chr / somme_wi1;
291 
292  muL2 = somme2_L / somme_wi2;
293  muCh2 = somme2_Chr / somme_wi2;
294 
295  // Compute standard deviation and mixed terms
296  for(unsigned int i = 0; i < tab_point1.size(); i++)
297  {
298  variance1_L += tab_wi1[i] * pow(tab_color1[i][0] - muL1, 2);
299  variance2_L += tab_wi2[i] * pow(tab_color2[i][0] - muL2, 2);
300 
301  dH_sq += tab_wi1[i] * tab_dH[i];
302  sL12 += tab_wi1[i] * (tab_color1[i][0] - muL1) *
303  (tab_color2[i][0] - muL2); // covariance
304  }
305 
306  variance1_L = variance1_L / somme_wi1;
307  sL1 = (variance1_L < 0) ? 0 : sqrt(variance1_L);
308 
309  variance2_L = variance2_L / somme_wi2;
310  sL2 = (variance2_L < 0) ? 0 : sqrt(variance2_L);
311 
312  sL12 = sL12 / somme_wi1;
313  sL12 = (sL12 < 0) ? 0 : sL12;
314 
315  dL_sq = (muL1 - muL2) * (muL1 - muL2);
316  dCh_sq = (muCh1 - muCh2) * (muCh1 - muCh2);
317 
318  dH_sq = dH_sq / somme_wi1;
319  dH_sq = dH_sq * dH_sq;
320 
321  // We then compute the local color-based features
322  double c1 = 0.002, c2 = 0.1, c3 = 0.1, c4 = 0.002,
323  c5 = 0.008; // The constants
324  // 1) Lightness comparison
325  double f4 = 1 - (1 / (c1 * dL_sq + 1));
326  // 2) Lightness-contrast comparison
327  double f5 = 1 - ((2 * sL1 * sL2 + c2) / (sL1 * sL1 + sL2 * sL2 + c2));
328  // 3) Lightness-structure comparison
329  double f6 = 1 - ((sL12 + c3) / (sL1 * sL2 + c3));
330  // 4) Chroma comparison
331  double f7 = 1 - (1 / (c4 * dCh_sq + 1));
332  // 5) Hue comparison
333  double f8 = 1 - (1 / (c5 * dH_sq + 1));
334 
335  tab_f4[p_vertex] += f4;
336  tab_f5[p_vertex] += f5;
337  tab_f6[p_vertex] += f6;
338  tab_f7[p_vertex] += f7;
339  tab_f8[p_vertex] += f8;
340 }
341 
351 template< typename HalfedgeGraph,
352  typename VertexColorMap = typename FEVV::
353  PMap_traits< FEVV::vertex_color_t, MeshSurface >::pmap_type >
354 VertexColorMap
356  const HalfedgeGraph &m_poly,
357  FEVV::PMapsContainer &pmaps,
358  std::vector< double > &init_grid_L,
359  std::vector< double > &grid_L,
360  std::vector< std::vector< std::pair< double, double > > > &init_grid_AB)
361 {
362  VertexColorMap v_cm;
363  v_cm = get_property_map(FEVV::vertex_color, m_poly, pmaps);
364 
365  // Transform the color of each vertex into the LAB2000HL space
366  // RGB -> LAB -> LAB2000HL
367  using Color = typename boost::property_traits< VertexColorMap >::value_type;
368  VertexColorMap _LAB2000HL_Vertexcolormap;
369 
370  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly))
371  {
372  double *RGB_Vcolor = new double[3]; // vertex color in RGB space
373  double *LAB_Vcolor = new double[3]; // vertex color in LAB space
374  double *LAB2000HL_Vcolor = new double[3]; // vertex color in LAB2000HL space
375 
376  RGB_Vcolor[0] = (v_cm)[vi][0];
377  RGB_Vcolor[1] = (v_cm)[vi][1];
378  RGB_Vcolor[2] = (v_cm)[vi][2];
379 
380  // Color space conversion.
381  rgb_to_lab(RGB_Vcolor[0], RGB_Vcolor[1], RGB_Vcolor[2], LAB_Vcolor);
382  lab2000hl_conversion(LAB_Vcolor[0],
383  LAB_Vcolor[1],
384  LAB_Vcolor[2],
385  LAB2000HL_Vcolor,
386  init_grid_L,
387  grid_L,
388  init_grid_AB);
389 
390  put(_LAB2000HL_Vertexcolormap,
391  vi,
392  Color(LAB2000HL_Vcolor[0], LAB2000HL_Vcolor[1], LAB2000HL_Vcolor[2]));
393 
394  delete[] RGB_Vcolor;
395  delete[] LAB_Vcolor;
396  delete[] LAB2000HL_Vcolor;
397  }
398 
399  return _LAB2000HL_Vertexcolormap;
400 }
401 
411 template< typename CMDMMap,
412  typename HalfedgeGraph,
413  typename CMDMNearestMap,
414  typename TagMap,
415  typename PointMap,
416  typename FaceIterator =
417  typename boost::graph_traits< HalfedgeGraph >::face_iterator >
418 void
420  const HalfedgeGraph &m_poly_degrad,
421  CMDMMap &cmdm_degrad,
422  CMDMNearestMap &cmdmm_degrad,
423  TagMap &tag_map,
424  const PointMap &pm_degrad, // Or CGALPoint *tab_pm_degrad,
425  const HalfedgeGraph &m_poly_original,
426  CGAL::SM_Face_index *tab_matched_facet)
427 {
428  using AABB_Primitive =
429  typename CGAL::AABB_face_graph_triangle_primitive< HalfedgeGraph >;
430  using AABB_Traits = typename CGAL::AABB_traits< CGALKernel, AABB_Primitive >;
431  using AABB_Tree = typename CGAL::AABB_tree< AABB_Traits >;
432  //using Object_and_primitive_id = typename AABB_Tree::Object_and_primitive_id;
433  using Point_and_primitive_id = typename AABB_Tree::Point_and_primitive_id;
434 
435  auto original_faces = faces(m_poly_original);
436 
437  FaceIterator orig_begin = original_faces.first;
438  FaceIterator orig_end = original_faces.second;
439 
440  AABB_Tree tree(orig_begin, orig_end, m_poly_original);
441  tree.accelerate_distance_queries();
442 
443  // Searching for the closest point and facet for each vertex
444 
445  // CGALPoint *tab_cmdmm_degrad = new CGALPoint[(int)num_vertices(m_poly_degrad)]; - used in parallelization trials
446 
447  //#pragma omp parallel for - program crashes with a critical error
448  // detected c0000374 (problem with Parallelism of
449  // AABB_tree::closest_point_and_primitive() function)
450  for(int i = 0; i < (int)num_vertices(m_poly_degrad); i++) // BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
451  {
452  Vertex_Descriptor vi(i);
453  // computes closest point and primitive id
454  Point_and_primitive_id pp =
455  tree.closest_point_and_primitive(get(pm_degrad, vi));
456  auto f_nearest = pp.second; // closest primitive id
457 
458  /* formulation used in parallelization trials
459  auto *p = &tab_pm_degrad[vi];
460  Point_and_primitive_id pp =(tree.closest_point_and_primitive(*p));
461  tab_cmdmm_degrad[i] = pp.first;
462  */
463 
464  put(cmdmm_degrad, vi, pp.first);
465  tab_matched_facet[i] = f_nearest;
466  }
467 
468  int ind = 0;
469  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
470  {
471  put(cmdm_degrad, vi, 0.);
472  put(tag_map, vi, ind);
473  // put(cmdmm_degrad, vi, tab_cmdmm_degrad[ind]);
474  ind++;
475  }
476 }
477 
478 
479 template< typename PointMap,
480  typename HalfedgeGraph,
481  typename FaceNormalMap,
482  typename VertexCMDMMap >
483 double
484 process_CMDM_multires(const HalfedgeGraph &m_poly_degrad,
485  const PointMap &pm_degrad,
486  const FaceNormalMap &fnm_degrad,
487  FEVV::PMapsContainer &pmaps_degrad,
488  const HalfedgeGraph &m_poly_original,
489  const PointMap &pm_original,
490  const FaceNormalMap &fnm_original,
491  FEVV::PMapsContainer &pmaps_original,
492  const int nb_level,
493  VertexCMDMMap &cmdm_pmap,
494  const double maxdim,
495  double &CMDM_value)
496 {
497  auto curvature_map_degr =
498  FEVV::make_vertex_property_map< HalfedgeGraph,
500  m_poly_degrad);
501  auto curvature_map_orig =
502  FEVV::make_vertex_property_map< HalfedgeGraph,
504  m_poly_original);
505 
506  CGAL::SM_Face_index *tab_matched_facet =
507  new CGAL::SM_Face_index[size_of_vertices(m_poly_degrad)];
508 
509 
510  VertexCMDMMap curv_match_pmap;
511 
512  if(has_map(pmaps_degrad, std::string("v:cmdm_match")))
513  {
514  curv_match_pmap =
515  boost::any_cast< VertexCMDMMap >(pmaps_degrad.at("v:cmdm_match"));
516  }
517  else
518  {
519  curv_match_pmap =
520  FEVV::make_vertex_property_map< HalfedgeGraph, double >(m_poly_degrad);
521  pmaps_degrad["v:cmd_match"] = curv_match_pmap;
522  }
523 
525  CMDM_nearest_map;
526  CMDM_nearest_map cmdm_nearest_pmap;
527 
528  if(has_map(pmaps_degrad, std::string("v:cmdm_nearest")))
529  {
530  cmdm_nearest_pmap = boost::any_cast< CMDM_nearest_map >(
531  pmaps_degrad.at("v:cmdm_nearest"));
532  }
533  else
534  {
535  cmdm_nearest_pmap =
536  FEVV::make_vertex_property_map< HalfedgeGraph, CGALPoint >(
537  m_poly_degrad);
538  pmaps_degrad["v:cmdm_nearest"] = cmdm_nearest_pmap;
539  }
540 
541  typedef typename FEVV::Vertex_pmap< HalfedgeGraph, int > vertex_tag_map;
542  vertex_tag_map tag_map;
543 
544  if(has_map(pmaps_degrad, std::string("v:tag")))
545  {
546  tag_map = boost::any_cast< vertex_tag_map >(pmaps_degrad.at("v:tag"));
547  }
548  else
549  {
550  tag_map =
551  FEVV::make_vertex_property_map< HalfedgeGraph, int >(m_poly_degrad);
552  pmaps_degrad["v:tag"] = tag_map;
553  }
554 
555  /***Note that, static arrays were used to store feature values and vertices
556  coordinates instead of proprty maps (PointMap, CMDMMap; e.g. CMDMMap
557  &feature1), as the latter caused problems when parallelizing the code (notably
558  the get() and put() functions): get(pm_degrad, vi) => tab_pm_degrad[i] ****/
559  double *tab_f1 = new double[(int)num_vertices(m_poly_degrad)](); // initialized to 0
560  double *tab_f2 = new double[(int)num_vertices(m_poly_degrad)]();
561  double *tab_f3 = new double[(int)num_vertices(m_poly_degrad)]();
562  double *tab_f4 = new double[(int)num_vertices(m_poly_degrad)]();
563  double *tab_f5 = new double[(int)num_vertices(m_poly_degrad)]();
564  double *tab_f6 = new double[(int)num_vertices(m_poly_degrad)]();
565  double *tab_f7 = new double[(int)num_vertices(m_poly_degrad)]();
566  double *tab_f8 = new double[(int)num_vertices(m_poly_degrad)]();
567 
568  CGALPoint *tab_pm_original =
569  new CGALPoint[(int)num_vertices(m_poly_original)];
570  CGALPoint *tab_pm_degrad =
571  new CGALPoint[(int)num_vertices(m_poly_degrad)];
572  CGALPoint *tab_cmdm_nearest_pmap =
573  new CGALPoint[(int)num_vertices(m_poly_degrad)];
574 
575 //#pragma omp parallel for
576  for(int ii = 0; ii < (int)num_vertices(m_poly_original); ii++)
577  {
578  Vertex_Descriptor vi(ii);
579  auto p0 = get(pm_original, vi);
580  tab_pm_original[ii] = p0;
581  }
582 
583 //#pragma omp parallel for
584  for(int ii = 0; ii < (int)num_vertices(m_poly_degrad); ii++)
585  {
586  Vertex_Descriptor vi(ii);
587  auto p1 = get(pm_degrad, vi);
588  tab_pm_degrad[ii] = p1;
589  }
590 
591  matching_multires_init_cmdm(m_poly_degrad,
592  cmdm_pmap,
593  cmdm_nearest_pmap,
594  tag_map,
595  pm_degrad,
596  m_poly_original,
597  tab_matched_facet);
598 
599 //#pragma omp parallel for
600  for(int ii = 0; ii < (int)num_vertices(m_poly_degrad); ii++)
601  {
602  Vertex_Descriptor vi(ii);
603  auto p2 = get(cmdm_nearest_pmap, vi);
604  tab_cmdm_nearest_pmap[ii] = p2;
605  }
606 
607  // Initialize, from files, the 1D and 2D grids for LAB to LAB2000HL
608  // conversion/interpolation
609  static std::vector< double > init_grid_L;
610  static std::vector< double > grid_L;
611  static std::vector< std::vector< std::pair< double, double > > > init_grid_AB;
612  initMatLABCH(init_grid_L, grid_L, init_grid_AB);
613 
614  using VertexColorMap = typename FEVV::PMap_traits< FEVV::vertex_color_t,
615  MeshSurface >::pmap_type;
616  VertexColorMap v_cm_orig = Get_VertexColorMap(
617  m_poly_original, pmaps_original, init_grid_L, grid_L, init_grid_AB);
618  VertexColorMap v_cm_degr = Get_VertexColorMap(
619  m_poly_degrad, pmaps_degrad, init_grid_L, grid_L, init_grid_AB);
620  VertexColorMap v_cm_match_pmap; // empty map which will be filled by the color
621  // of the projected/matched points, in
622  // matching_multires_update() function
623 
624  double radius_curvature = 0.001;
625  for(int l = 0; l < nb_level; l++)
626  {
627  double min_nrm_min_curvature_deg, max_nrm_min_curvature_deg,
628  min_nrm_max_curvature_deg, max_nrm_max_curvature_deg;
629  double min_nrm_min_curvature_org, max_nrm_min_curvature_org,
630  min_nrm_max_curvature_org, max_nrm_max_curvature_org;
631 
633  curvature_map_degr,
634  pm_degrad,
635  fnm_degrad,
636  true,
637  maxdim * radius_curvature,
638  min_nrm_min_curvature_deg,
639  max_nrm_min_curvature_deg,
640  min_nrm_max_curvature_deg,
641  max_nrm_max_curvature_deg);
643  curvature_map_orig,
644  pm_original,
645  fnm_original,
646  true,
647  maxdim * radius_curvature,
648  min_nrm_min_curvature_org,
649  max_nrm_min_curvature_org,
650  min_nrm_max_curvature_org,
651  max_nrm_max_curvature_org);
652 
653  kmax_kmean_cmdm(m_poly_original, curvature_map_orig, maxdim);
654  kmax_kmean_cmdm(m_poly_degrad, curvature_map_degr, maxdim);
655 
656  matching_multires_update(m_poly_degrad,
657  m_poly_original,
658  tab_pm_original,
659  tab_cmdm_nearest_pmap,
660  curvature_map_orig,
661  v_cm_orig,
662  curv_match_pmap,
663  v_cm_match_pmap,
664  tab_matched_facet);
665 
666  //#pragma omp parallel for
667  for(int i = 0; i < (int)num_vertices(m_poly_degrad);
668  ++i) // for(Vertex_Descriptor vi: vertices(m_poly_degrad))
669  {
670  // std::cout << "nb thread " << omp_get_num_threads() << std::endl;
671  Vertex_Descriptor vi(i);
672 
673  std::vector< CGALPoint > tab_point1;
674  std::vector< CGALPoint > tab_point2;
675 
676  std::vector< double > tab_distance1;
677  std::vector< double > tab_distance2;
678 
679  std::vector< std::array< double, 3 > > tab_color1;
680  std::vector< std::array< double, 3 > > tab_color2;
681 
682  process_cmdm_per_vertex(m_poly_degrad,
683  tab_pm_degrad,
684  tag_map,
685  tab_cmdm_nearest_pmap,
686  curvature_map_degr,
687  v_cm_degr,
688  curv_match_pmap,
689  v_cm_match_pmap,
690  vi,
691  radius_curvature * 3 * maxdim,
692  tab_point1,
693  tab_point2,
694  tab_distance1,
695  tab_distance2,
696  tab_color1,
697  tab_color2);
698 
699  compute_geometry_statistics(m_poly_degrad,
700  pm_degrad,
701  cmdm_nearest_pmap,
702  vi,
703  tab_distance1,
704  tab_distance2,
705  tab_point1,
706  tab_point2,
707  tab_f1,
708  tab_f2,
709  tab_f3,
710  radius_curvature * 3 * maxdim);
711 
712  compute_color_statistics(m_poly_degrad,
713  pm_degrad,
714  cmdm_nearest_pmap,
715  vi,
716  tab_color1,
717  tab_color2,
718  tab_point1,
719  tab_point2,
720  tab_f4,
721  tab_f5,
722  tab_f6,
723  tab_f7,
724  tab_f8,
725  radius_curvature * 3 * maxdim);
726  }
727  radius_curvature += 0.0005;
728  }
729 
730  double sum_f1 = 0, sum_f2 = 0, sum_f3 = 0, sum_f4 = 0, sum_f5 = 0, sum_f6 = 0,
731  sum_f7 = 0, sum_f8 = 0;
732  double w2 = 0.091, w5 = 0.22, w6 = 0.032,
733  w7 = 0.656; // The recommended weights
734  int nb_vertex = (int)num_vertices(m_poly_degrad);
735 
736  //#pragma omp parallel for
737  for(int i = 0; i < nb_vertex; ++i)
738  {
739  //#pragma omp atomic
740  sum_f1 += tab_f1[i];
741  sum_f2 += tab_f2[i];
742  sum_f3 += tab_f3[i];
743  sum_f4 += tab_f4[i];
744  sum_f5 += tab_f5[i];
745  sum_f6 += tab_f6[i];
746  sum_f7 += tab_f7[i];
747  sum_f8 += tab_f8[i];
748 
749  // Local distortion
750  double LD =
751  (w2 * tab_f2[i] + w5 * tab_f5[i] + w6 * tab_f6[i] + w7 * tab_f7[i]) /
752  nb_level;
753  Vertex_Descriptor vi(i);
754  put(cmdm_pmap, vi, LD);
755  }
756 
757  sum_f1 = sum_f1 / (nb_vertex * nb_level);
758  sum_f2 = sum_f2 / (nb_vertex * nb_level);
759  sum_f3 = sum_f3 / (nb_vertex * nb_level);
760  sum_f4 = sum_f4 / (nb_vertex * nb_level);
761  sum_f5 = sum_f5 / (nb_vertex * nb_level);
762  sum_f6 = sum_f6 / (nb_vertex * nb_level);
763  sum_f7 = sum_f7 / (nb_vertex * nb_level);
764  sum_f8 = sum_f8 / (nb_vertex * nb_level);
765 
766  CMDM_value = w2 * sum_f2 + w5 * sum_f5 + w6 * sum_f6 + w7 * sum_f7;
767 
768  delete[] tab_pm_original;
769  delete[] tab_pm_degrad;
770  delete[] tab_matched_facet;
771  delete[] tab_cmdm_nearest_pmap;
772  delete[] tab_f1;
773  delete[] tab_f2;
774  delete[] tab_f3;
775  delete[] tab_f4;
776  delete[] tab_f5;
777  delete[] tab_f6;
778  delete[] tab_f7;
779  delete[] tab_f8;
780 
781  return 0;
782 }
783 
792 template< typename HalfedgeGraph,
793  typename CurvatureVertexMap,
794  typename CMDMCurvatureMatchMap,
795  typename VertexColorMap >
796 void
798  const HalfedgeGraph &m_poly_degrad,
799  const HalfedgeGraph &m_poly_orig,
800  CGALPoint *tab_pm_orig, // Or const PointMap &pm_orig,
801  CGALPoint *tab_cmdmm_degrad, // Or CMDMNearestMap &cmdmm_degrad,
802  const CurvatureVertexMap curv_orig,
803  const VertexColorMap vcm_orig,
804  CMDMCurvatureMatchMap curvmatch_pmap,
805  VertexColorMap colmatch_pmap,
806  CGAL::SM_Face_index *tab_matched_facet)
807 {
808  using Color = typename boost::property_traits< VertexColorMap >::value_type;
809  Color *tab_colmatch_pmap = new Color[(int)num_vertices(m_poly_degrad)];
810  double *tab_curvmatch_pmap = new double[(int)num_vertices(m_poly_degrad)];
811 
812  //#pragma omp parallel for
813  for(int i = 0; i < (int)num_vertices(m_poly_degrad);
814  ++i) // BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
815  {
816  // Vertex_Descriptor vi(i);
817  int ind = i;
818  auto *f_nearest = &tab_matched_facet[ind];
819 
822  //we use linear interpolation using barycentric coordinates
823 
824  /* formulation used in NOT parallelized function
825  auto x1 = get(pm_orig, target(halfedge(*f_nearest, m_poly_orig),
826  m_poly_orig)); auto x2 = get(pm_orig, target(next(halfedge(*f_nearest,
827  m_poly_orig), m_poly_orig), m_poly_orig)); auto x3 =get(pm_orig,
828  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
829  m_poly_orig), m_poly_orig));
830  */
831 
832  auto *x1 =
833  &tab_pm_orig[target(halfedge(*f_nearest, m_poly_orig), m_poly_orig)];
834  auto *x2 = &tab_pm_orig[target(
835  next(halfedge(*f_nearest, m_poly_orig), m_poly_orig), m_poly_orig)];
836  auto *x3 = &tab_pm_orig[target(
837  next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig), m_poly_orig),
838  m_poly_orig)];
839 
840  double l1 = sqrt((*x3 - *x2) * (*x3 - *x2));
841  double l2 = sqrt((*x1 - *x3) * (*x1 - *x3));
842  double l3 = sqrt((*x1 - *x2) * (*x1 - *x2));
843 
844  auto *pt_matched = &tab_cmdmm_degrad[i];
845  auto v1 = *x1 - *pt_matched; // x1 - get(cmdmm_degrad, vi);
846  auto v2 = *x2 - *pt_matched; // x2 - get(cmdmm_degrad, vi);
847  auto v3 = *x3 - *pt_matched; // x3 - get(cmdmm_degrad, vi);
848 
849  double t1 = sqrt(v1 * v1);
850  double t2 = sqrt(v2 * v2);
851  double t3 = sqrt(v3 * v3);
852 
853  double p1 = (l1 + t2 + t3) / 2;
854  double p2 = (t1 + l2 + t3) / 2;
855  double p3 = (t1 + t2 + l3) / 2;
856 
857  double a1 = (p1 * (p1 - l1) * (p1 - t3) * (p1 - t2));
858  double a2 = (p2 * (p2 - l2) * (p2 - t3) * (p2 - t1));
859  double a3 = (p3 * (p3 - l3) * (p3 - t1) * (p3 - t2));
860 
861  if(a1 > 0)
862  a1 = sqrt(a1);
863  else
864  a1 = 0;
865  if(a2 > 0)
866  a2 = sqrt(a2);
867  else
868  a2 = 0;
869  if(a3 > 0)
870  a3 = sqrt(a3);
871  else
872  a3 = 0;
873 
874  double c1 =
875  get(curv_orig, target(halfedge(*f_nearest, m_poly_orig), m_poly_orig))
876  .KmaxCurv;
877  double c2 = get(curv_orig,
878  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
879  m_poly_orig))
880  .KmaxCurv;
881  double c3 =
882  get(curv_orig,
883  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
884  m_poly_orig),
885  m_poly_orig))
886  .KmaxCurv;
887 
888  std::array< double, 3 > col1 = {
889  get(vcm_orig,
890  target(halfedge(*f_nearest, m_poly_orig), m_poly_orig))[0],
891  get(vcm_orig,
892  target(halfedge(*f_nearest, m_poly_orig), m_poly_orig))[1],
893  get(vcm_orig,
894  target(halfedge(*f_nearest, m_poly_orig), m_poly_orig))[2]};
895  std::array< double, 3 > col2 = {
896  get(vcm_orig,
897  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
898  m_poly_orig))[0],
899  get(vcm_orig,
900  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
901  m_poly_orig))[1],
902  get(vcm_orig,
903  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
904  m_poly_orig))[2]};
905  std::array< double, 3 > col3 = {
906  get(vcm_orig,
907  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
908  m_poly_orig),
909  m_poly_orig))[0],
910  get(vcm_orig,
911  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
912  m_poly_orig),
913  m_poly_orig))[1],
914  get(vcm_orig,
915  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
916  m_poly_orig),
917  m_poly_orig))[2]};
918 
919  double curvmatch = 0.0;
920  std::array< double, 3 > colmatch = {0, 0, 0};
921  if((a1 + a2 + a3) > 0)
922  {
923  curvmatch = (a1 * c1 + a2 * c2 + a3 * c3) / (a1 + a2 + a3);
924  colmatch = {(a1 * col1[0] + a2 * col2[0] + a3 * col3[0]) / (a1 + a2 + a3),
925  (a1 * col1[1] + a2 * col2[1] + a3 * col3[1]) / (a1 + a2 + a3),
926  (a1 * col1[2] + a2 * col2[2] + a3 * col3[2]) /
927  (a1 + a2 + a3)};
928  }
929  else
930  {
931  curvmatch = (c1 + c2 + c3) / 3;
932  colmatch = {(col1[0] + col2[0] + col3[0]) / 3,
933  (col1[1] + col2[1] + col3[1]) / 3,
934  (col1[2] + col2[2] + col3[2]) / 3};
935  }
936 
937  tab_curvmatch_pmap[i] = curvmatch; // put(curvmatch_pmap, vi, curvmatch);
938  tab_colmatch_pmap[i] =
939  Color(colmatch[0],
940  colmatch[1],
941  colmatch[2]); // put(colmatch_pmap, vi, Color(colmatch[0],
942  // colmatch[1], colmatch[2]));
943  }
944 
945  int index = 0;
946  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
947  {
948  put(curvmatch_pmap, vi, tab_curvmatch_pmap[index]);
949  put(colmatch_pmap, vi, tab_colmatch_pmap[index]);
950  index++;
951  }
952 
953  delete[] tab_curvmatch_pmap;
954  delete[] tab_colmatch_pmap;
955 }
956 
957 
977 template< typename VertexIterator,
978  typename HalfedgeGraph,
979  typename TagMap,
980  typename CurvatureVertexMap,
981  typename CMDMCurvatureMatchMap,
982  typename VertexColorMap >
983 void
985  const HalfedgeGraph &m_degr,
986  CGALPoint *tab_pm_degr, // const PointMap &pm_degr,
987  const TagMap &tags_pmap,
988  CGALPoint *tab_nearest_pmap, // const CMDMNearestMap &nearest_pmap,
989  const CurvatureVertexMap &curv_pmap_degr,
990  const VertexColorMap &col_pmap_degr,
991  const CMDMCurvatureMatchMap &curvmatch_pmap,
992  const VertexColorMap &colmatch_pmap,
993  const VertexIterator &p_vertex,
994  double radius,
995  std::vector< CGALPoint > &tab_point1,
996  std::vector< CGALPoint > &tab_point2,
997  std::vector< double > &tab_distance1,
998  std::vector< double > &tab_distance2,
999  std::vector< std::array< double, 3 > > &tab_color1,
1000  std::vector< std::array< double, 3 > > &tab_color2)
1001 {
1002  std::set< int > vertices;
1003 
1004  std::stack< VertexIterator > s;
1005 
1006  auto *o = &tab_pm_degr[p_vertex]; // auto o = get(pm_degr, p_vertex);
1007 
1008  s.push(p_vertex);
1009  vertices.insert(get(tags_pmap, p_vertex));
1010 
1011 
1012  tab_point1.push_back(
1013  tab_pm_degr[p_vertex]); // tab_point1.push_back(get(pm_degr, p_vertex));
1014  tab_distance1.push_back(get(curv_pmap_degr, p_vertex).KmaxCurv);
1015  std::array< double, 3 > VertexColor1 = {(col_pmap_degr)[p_vertex][0],
1016  (col_pmap_degr)[p_vertex][1],
1017  (col_pmap_degr)[p_vertex][2]};
1018  tab_color1.push_back(VertexColor1);
1019 
1020  tab_point2.push_back(tab_nearest_pmap[p_vertex]); // tab_point2.push_back(get(nearest_pmap,
1021  // p_vertex));
1022  tab_distance2.push_back(get(curvmatch_pmap, p_vertex));
1023  std::array< double, 3 > VertexColor2 = {(colmatch_pmap)[p_vertex][0],
1024  (colmatch_pmap)[p_vertex][1],
1025  (colmatch_pmap)[p_vertex][2]};
1026  tab_color2.push_back(VertexColor2);
1027 
1028  int nb_sommet_in_sphere = 0;
1029 
1030  while(!s.empty())
1031  {
1032  VertexIterator v_it = s.top();
1033  s.pop();
1034  CGALPoint *p = &tab_pm_degr[v_it]; // CGALPoint p = get(pm_degr, v_it);
1035  auto h = halfedge(v_it, m_degr);
1036  auto h0 = h;
1037  do
1038  {
1039  // A halfedge is directed from its source vertex to its target vertex
1040  auto v_target_of_h = target(h, m_degr);
1041  auto h_opp = opposite(h, m_degr);
1042  auto v_target_of_hopp = target(h_opp, m_degr);
1043 
1044  auto *p1 = &tab_pm_degr[v_target_of_h]; // auto p1 = get(pm_degr,
1045  // target(h, m_degr));
1046  auto *p2 = &tab_pm_degr[v_target_of_hopp]; // auto p2 = get(pm_degr,
1047  // target(opposite(h, m_degr),
1048  // m_degr));
1049 
1050  auto *p1m =
1051  &tab_nearest_pmap[v_target_of_h]; // auto p1m = get(nearest_pmap,
1052  // target(h, m_degr));
1053  auto *p2m =
1054  &tab_nearest_pmap[v_target_of_hopp]; // auto p2m = get(nearest_pmap,
1055  // target(opposite(h, m_degr),
1056  // m_degr));
1057 
1058  auto v = (*p2 - *p1);
1059  auto vm = (*p2m - *p1m);
1060 
1061  if(v_it == p_vertex || v * (*p - *o) >= 0.0)
1062  {
1063  double len_old = std::sqrt(v * v);
1064  bool isect = sphere_clip_vector_cmdm(*o, radius, *p, v);
1065  double len_edge = std::sqrt(v * v);
1066 
1067  nb_sommet_in_sphere++;
1068 
1069  CGALPoint weighted_p1, weighted_p2;
1070  double weighted_curv1, weighted_curv2;
1071  std::array< double, 3 > weighted_col1, weighted_col2;
1072 
1073  bool is_already_integrated = false;
1074  if(!isect)
1075  {
1076 
1077  auto w = target(opposite(h, m_degr), m_degr);
1078  if(vertices.find(get(tags_pmap, w)) == vertices.end())
1079  {
1080  vertices.insert(get(tags_pmap, w));
1081  s.push(w);
1082  }
1083  else
1084  is_already_integrated = true;
1085  }
1086 
1087  if(is_already_integrated == false)
1088  {
1089  if(len_old != 0 && isect)
1090  {
1091  weighted_p1 = *p1 + v;
1092  weighted_curv1 =
1093  (1 - len_edge / len_old) *
1094  get(curv_pmap_degr, target(h, m_degr)).KmaxCurv +
1095  len_edge / len_old *
1096  get(curv_pmap_degr, target(opposite(h, m_degr), m_degr))
1097  .KmaxCurv;
1098  weighted_col1 = {
1099  (1 - len_edge / len_old) *
1100  get(col_pmap_degr, target(h, m_degr))[0] +
1101  len_edge / len_old *
1102  get(col_pmap_degr,
1103  target(opposite(h, m_degr), m_degr))[0],
1104  (1 - len_edge / len_old) *
1105  get(col_pmap_degr, target(h, m_degr))[1] +
1106  len_edge / len_old *
1107  get(col_pmap_degr,
1108  target(opposite(h, m_degr), m_degr))[1],
1109  (1 - len_edge / len_old) *
1110  get(col_pmap_degr, target(h, m_degr))[2] +
1111  len_edge / len_old *
1112  get(col_pmap_degr,
1113  target(opposite(h, m_degr), m_degr))[2]};
1114 
1115  weighted_p2 = *p1m + (len_edge / len_old) * vm;
1116  weighted_curv2 =
1117  (1 - len_edge / len_old) *
1118  get(curvmatch_pmap, target(h, m_degr)) +
1119  len_edge / len_old *
1120  get(curvmatch_pmap, target(opposite(h, m_degr), m_degr));
1121  weighted_col2 = {
1122  (1 - len_edge / len_old) *
1123  get(colmatch_pmap, target(h, m_degr))[0] +
1124  len_edge / len_old *
1125  get(colmatch_pmap,
1126  target(opposite(h, m_degr), m_degr))[0],
1127  (1 - len_edge / len_old) *
1128  get(colmatch_pmap, target(h, m_degr))[1] +
1129  len_edge / len_old *
1130  get(colmatch_pmap,
1131  target(opposite(h, m_degr), m_degr))[1],
1132  (1 - len_edge / len_old) *
1133  get(colmatch_pmap, target(h, m_degr))[2] +
1134  len_edge / len_old *
1135  get(colmatch_pmap,
1136  target(opposite(h, m_degr), m_degr))[2]};
1137  }
1138  else
1139  {
1140  weighted_p1 = *p2;
1141  weighted_p2 = *p2m;
1142  weighted_curv1 =
1143  get(curv_pmap_degr, target(opposite(h, m_degr), m_degr))
1144  .KmaxCurv;
1145  weighted_curv2 =
1146  get(curvmatch_pmap, target(opposite(h, m_degr), m_degr));
1147  weighted_col1 = {
1148  get(col_pmap_degr, target(opposite(h, m_degr), m_degr))[0],
1149  get(col_pmap_degr, target(opposite(h, m_degr), m_degr))[1],
1150  get(col_pmap_degr, target(opposite(h, m_degr), m_degr))[2]};
1151  weighted_col2 = {
1152  get(colmatch_pmap, target(opposite(h, m_degr), m_degr))[0],
1153  get(colmatch_pmap, target(opposite(h, m_degr), m_degr))[1],
1154  get(colmatch_pmap, target(opposite(h, m_degr), m_degr))[2]};
1155  }
1156 
1157  tab_point1.push_back(weighted_p1);
1158  tab_distance1.push_back(weighted_curv1);
1159  tab_color1.push_back(weighted_col1);
1160 
1161  tab_point2.push_back(weighted_p2);
1162  tab_distance2.push_back(weighted_curv2);
1163  tab_color2.push_back(weighted_col2);
1164  }
1165  }
1166  h = opposite(next(h, m_degr), m_degr);
1167  } while(h != h0);
1168  }
1169 }
1170 
1171 
1179 template< typename HalfedgeGraph, typename CurvatureVertexMap >
1180 void
1181 kmax_kmean_cmdm(const HalfedgeGraph &mesh,
1182  CurvatureVertexMap &curv_vmap,
1183  const double coef)
1184 {
1185  BOOST_FOREACH(Vertex_Descriptor vi, vertices(mesh))
1186  {
1187  auto curv = get(curv_vmap, vi);
1188  double kmax = curv.KmaxCurv * coef;
1189  double kmin = std::abs(curv.KminCurv * coef);
1190  curv.KmaxCurv = (kmax + kmin) / 2.;
1191  put(curv_vmap, vi, curv);
1192  }
1193 }
1194 
1195 
1205 template< typename HalfedgeGraph,
1206  typename GeometryTraits = FEVV::Geometry_traits< HalfedgeGraph >,
1207  typename PointMap,
1208  typename FaceNormalMap,
1209  typename VertexCMDMMap >
1210 double
1211 process_CMDM_multires(const HalfedgeGraph &m_poly_degrad,
1212  const PointMap &pm_degrad,
1213  const FaceNormalMap &fnm_degrad,
1214  FEVV::PMapsContainer &pmaps_degrad,
1215  const HalfedgeGraph &m_poly_original,
1216  const PointMap &pm_original,
1217  const FaceNormalMap &fnm_original,
1218  FEVV::PMapsContainer &pmaps_original,
1219  const int nb_level,
1220  VertexCMDMMap &cmdm_pmap,
1221  double &CMDM_value)
1222 {
1223  GeometryTraits gt_deg(m_poly_degrad);
1224  GeometryTraits gt_org(m_poly_original);
1225 
1226  return process_CMDM_multires(
1227  m_poly_degrad,
1228  pm_degrad,
1229  fnm_degrad,
1230  pmaps_degrad,
1231  m_poly_original,
1232  pm_original,
1233  fnm_original,
1234  pmaps_original,
1235  nb_level,
1236  cmdm_pmap,
1237  Filters::get_max_bb_size(m_poly_degrad, pm_degrad, gt_deg),
1238  CMDM_value);
1239 }
1240 
1241 } // namespace Filters
1242 } // namespace FEVV
FEVV::DataStructures::AIF::vertices
std::pair< typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::vertex_iterator, typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::vertex_iterator > vertices(const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns the iterator range of the vertices of the mesh.
Definition: Graph_traits_aif.h:172
FEVV::Math::Vector::Stats::variance
static ElementType variance(const std::vector< ElementType > &v, ElementType mean_value, bool unbiased=true)
Definition: MatrixOperations.hpp:354
FEVV::DataStructures::AIF::num_vertices
boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::vertices_size_type num_vertices(const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns an upper bound of the number of vertices of the mesh.
Definition: Graph_traits_aif.h:191
FEVV::DataStructures::AIF::next
boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::halfedge_descriptor next(typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::halfedge_descriptor h, const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns the next halfedge around its face.
Definition: Graph_traits_aif.h:599
FEVV::get_property_map
PMap_traits< PropertyT, MeshT >::pmap_type get_property_map(PropertyT p, const MeshT &, const PMapsContainer &pmaps)
Definition: properties.h:646
FEVV::Filters::compute_geometry_statistics
void compute_geometry_statistics(const HalfedgeGraph &, const PointMap &pm_degr, const CMDMNearestMap &nearest_pmap, const VertexIterator p_vertex, std::vector< double > &tab_distance1, std::vector< double > &tab_distance2, std::vector< CGALPoint > &tab_point1, std::vector< CGALPoint > &tab_point2, double *tab_f1, double *tab_f2, double *tab_f3, double radius)
Calculates the local geometry-based features per vertex.
Definition: cmdm.h:112
rgb_to_lab.h
lab2000hl_conversion
void lab2000hl_conversion(double L, double A, double B, double *lab2000hl, std::vector< double > &init_grid_L, std::vector< double > &grid_L, std::vector< std::vector< std::pair< double, double > > > &init_grid_AB)
Definition: transfo_lab2lab2000hl.h:175
FEVV::Vertex_Descriptor
MeshSurface::Vertex_index Vertex_Descriptor
Definition: cmdm.h:44
FEVV::has_map
bool has_map(const PMapsContainer &pmaps, const std::string &map_name)
(refer to Property Maps API)
Definition: properties.h:103
generic_reader.hpp
FEVV::Filters::calculate_curvature_cmdm
void calculate_curvature_cmdm(const HalfedgeGraph &g, VertexCurvatureMap &v_cm, const PointMap &pm, const FaceNormalMap &f_nm, bool is_geod, double radius, double &min_nrm_min_curvature, double &max_nrm_min_curvature, double &min_nrm_max_curvature, double &max_nrm_max_curvature, const GeometryTraits &gt)
Calculate the curvature for a mesh.
Definition: curvature_cmdm.hpp:523
FEVV::DataStructures::AIF::opposite
boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::halfedge_descriptor opposite(typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::halfedge_descriptor h, const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns the halfedge with source and target swapped.
Definition: Graph_traits_aif.h:625
FEVV::Geometry_traits< HalfedgeGraph >
Geometry_traits_cgal_surface_mesh.h
initMatLABCH
void initMatLABCH(std::vector< double > &init_grid_L, std::vector< double > &grid_L, std::vector< std::vector< std::pair< double, double > > > &init_grid_AB)
Definition: transfo_lab2lab2000hl.h:137
FEVV::MeshSurface
CGAL::Surface_mesh< CGALPoint > MeshSurface
Definition: DataStructures_cgal_surface_mesh.h:23
FEVV::PMapsContainer
std::map< std::string, boost::any > PMapsContainer
Definition: properties.h:99
FEVV::Filters::matching_multires_init_cmdm
void matching_multires_init_cmdm(const HalfedgeGraph &m_poly_degrad, CMDMMap &cmdm_degrad, CMDMNearestMap &cmdmm_degrad, TagMap &tag_map, const PointMap &pm_degrad, const HalfedgeGraph &m_poly_original, CGAL::SM_Face_index *tab_matched_facet)
Initialize the matching process.
Definition: cmdm.h:419
FEVV::get
FEVV::PCLPointCloudPointMap::value_type get(const FEVV::PCLPointCloudPointMap &pm, FEVV::PCLPointCloudPointMap::key_type key)
Specialization of get(point_map, key) for PCLPointCloud.
Definition: Graph_properties_pcl_point_cloud.h:117
FEVV::Filters::kmax_kmean_cmdm
void kmax_kmean_cmdm(const HalfedgeGraph &mesh, CurvatureVertexMap &curv_vmap, const double coef)
Computes the mean curvature field (kmin+kmax)/2 and normalize it according to the size of the model.
Definition: cmdm.h:1181
FEVV::vertex_color_t
vertex_color_t
Definition: properties.h:47
FEVV
Interfaces for plugins These interfaces will be used for different plugins.
Definition: Assert.h:16
FEVV::DataStructures::AIF::faces
std::pair< typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::face_iterator, typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::face_iterator > faces(const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns an iterator range over all faces of the mesh.
Definition: Graph_traits_aif.h:679
FEVV::Filters::matching_multires_update
void matching_multires_update(const HalfedgeGraph &m_poly_degrad, const HalfedgeGraph &m_poly_orig, CGALPoint *tab_pm_orig, CGALPoint *tab_cmdmm_degrad, const CurvatureVertexMap curv_orig, const VertexColorMap vcm_orig, CMDMCurvatureMatchMap curvmatch_pmap, VertexColorMap colmatch_pmap, CGAL::SM_Face_index *tab_matched_facet)
Updates the matching process.
Definition: cmdm.h:797
properties_surface_mesh.h
FEVV::Filters::get_max_bb_size
double get_max_bb_size(const HalfedgeGraph &g, PointMap &pm, GeometryTraits &gt)
Gets the maximum size of the object bounding box.
Definition: get_max_bb_size.hpp:36
FEVV::DataStructures::AIF::halfedge
boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::halfedge_descriptor halfedge(typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::vertex_descriptor v, const FEVV::DataStructures::AIF::AIFMesh &sm)
Returns a halfedge with target v.
Definition: Graph_traits_aif.h:296
FEVV::DataStructures::AIF::target
boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::vertex_descriptor target(typename boost::graph_traits< FEVV::DataStructures::AIF::AIFMesh >::edge_descriptor e, const FEVV::DataStructures::AIF::AIFMesh &)
Returns the target vertex of e.
Definition: Graph_traits_aif.h:400
FEVV::Color
Definition: Color.hpp:18
AABB_Primitive
CGAL::AABB_triangle_primitive< AABB_Kernel, std::list< Triangle >::iterator > AABB_Primitive
A primitive for an AABB-tree.
Definition: boolops_polyhedra.hpp:121
FEVV::vertex_color
@ vertex_color
Definition: properties.h:47
FEVV::Ray
CGALKernel::Ray_3 Ray
Definition: cmdm.h:43
FEVV::make_vertex_property_map
Vertex_pmap_traits< MeshT, ValueT >::pmap_type make_vertex_property_map(const MeshT &m)
Definition: properties.h:736
FEVV::Filters::v_Curv_cmdm
Definition: curvature_cmdm.hpp:41
FEVV::DataStructures::AIF::AIFVector
Definition: AIFProperties.h:173
DataStructures_cgal_surface_mesh.h
FEVV::CGALPoint
CGALKernel::Point_3 CGALPoint
Definition: DataStructures_cgal_surface_mesh.h:21
FEVV::Vertex_pmap
typename Vertex_pmap_traits< MeshT, ValueT >::pmap_type Vertex_pmap
Definition: properties.h:607
FEVV::Filters::sphere_clip_vector_cmdm
bool sphere_clip_vector_cmdm(Point3d &o, double r, const Point3d &p, Vector &v)
Definition: cmdm.h:50
AABB_Traits
CGAL::AABB_traits< AABB_Kernel, AABB_Primitive > AABB_Traits
concept for AABB-tree
Definition: boolops_polyhedra.hpp:124
generic_writer.hpp
FEVV::Filters::fabs
double fabs(const v_Curv< HalfedgeGraph > &input)
Definition: curvature.hpp:54
FEVV::Triangle
CGALKernel::Triangle_3 Triangle
Definition: cmdm.h:41
FEVV::put
void put(FEVV::PCLPointCloudPointMap &pm, FEVV::PCLPointCloudPointMap::key_type key, const FEVV::PCLPointCloudPointMap::value_type &value)
Specialization of put(point_map, key, value) for PCLPointCloud.
Definition: Graph_properties_pcl_point_cloud.h:126
AABB_Tree
CGAL::AABB_tree< AABB_Traits > AABB_Tree
AABB-tree.
Definition: boolops_polyhedra.hpp:127
FEVV::Filters::process_cmdm_per_vertex
void process_cmdm_per_vertex(const HalfedgeGraph &m_degr, CGALPoint *tab_pm_degr, const TagMap &tags_pmap, CGALPoint *tab_nearest_pmap, const CurvatureVertexMap &curv_pmap_degr, const VertexColorMap &col_pmap_degr, const CMDMCurvatureMatchMap &curvmatch_pmap, const VertexColorMap &colmatch_pmap, const VertexIterator &p_vertex, double radius, std::vector< CGALPoint > &tab_point1, std::vector< CGALPoint > &tab_point2, std::vector< double > &tab_distance1, std::vector< double > &tab_distance2, std::vector< std::array< double, 3 > > &tab_color1, std::vector< std::array< double, 3 > > &tab_color2)
Computes the local neighborhoods.
Definition: cmdm.h:984
Point3d
EnrichedPolyhedron::Point_3 Point3d
Definition: boolops_cpolyhedron_builder.hpp:27
transfo_lab2lab2000hl.h
get_max_bb_size.hpp
FEVV::Segment
CGALKernel::Segment_3 Segment
Definition: cmdm.h:42
FEVV::Filters::Get_VertexColorMap
VertexColorMap Get_VertexColorMap(const HalfedgeGraph &m_poly, FEVV::PMapsContainer &pmaps, std::vector< double > &init_grid_L, std::vector< double > &grid_L, std::vector< std::vector< std::pair< double, double > > > &init_grid_AB)
Retrieve the color (RGB values) of each vertex of the mesh and transform it into the working color sp...
Definition: cmdm.h:355
FEVV::Filters::compute_color_statistics
void compute_color_statistics(const HalfedgeGraph &, const PointMap &pm_degr, const CMDMNearestMap &nearest_pmap, const VertexIterator p_vertex, std::vector< std::array< double, 3 > > &tab_color1, std::vector< std::array< double, 3 > > &tab_color2, std::vector< CGALPoint > &tab_point1, std::vector< CGALPoint > &tab_point2, double *tab_f4, double *tab_f5, double *tab_f6, double *tab_f7, double *tab_f8, double radius)
Calculates the local color-based features per vertex.
Definition: cmdm.h:214
rgb_to_lab
void rgb_to_lab(double R, double G, double B, double *lab)
Definition: rgb_to_lab.h:67
FEVV::_PMap_traits
Definition: properties.h:376
FEVV::Filters::process_CMDM_multires
double process_CMDM_multires(const HalfedgeGraph &m_poly_degrad, const PointMap &pm_degrad, const FaceNormalMap &fnm_degrad, FEVV::PMapsContainer &pmaps_degrad, const HalfedgeGraph &m_poly_original, const PointMap &pm_original, const FaceNormalMap &fnm_original, FEVV::PMapsContainer &pmaps_original, const int nb_level, VertexCMDMMap &cmdm_pmap, const double maxdim, double &CMDM_value)
Definition: cmdm.h:484
curvature_cmdm.hpp
FEVV::size_of_vertices
boost::graph_traits< MeshT >::vertices_size_type size_of_vertices(const MeshT &g)
Real current number of vertices of the mesh. Generic version.
Definition: Graph_traits_extension.h:29