MEPP2 Project
msdm2.h
Go to the documentation of this file.
1 // Copyright (c) 2012-2019 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 #pragma once
12 
14 
17 
18 #include <CGAL/AABB_tree.h>
19 #include <CGAL/AABB_traits.h>
20 #include <CGAL/AABB_face_graph_triangle_primitive.h>
21 
22 //---------------- specific code above --------------------
23 
24 //---------------- generic code below ---------------------
25 
26 // The code below is generic and works with all mesh datastructures
27 // without any change.
28 
33 
34 #include <boost/foreach.hpp>
35 
36 namespace FEVV {
37 
38 using Triangle = CGALKernel::Triangle_3;
39 using Segment = CGALKernel::Segment_3;
40 using Ray = CGALKernel::Ray_3;
41 using Vertex_Descriptor = MeshSurface::Vertex_index;
42 
43 namespace Filters {
44 
45 template< typename Point3d, typename Vector >
46 bool
47 sphere_clip_vector_msdm2(Point3d &o, double r, const Point3d &p, Vector &v)
48 {
49 
50  Vector w = p - o;
51  double a = (v * v);
52  double b = 2.0 * v * w;
53  double c = (w * w) - r * r;
54  double delta = b * b - 4 * a * c;
55 
56  if(a == 0)
57  return true;
58 
59  if(delta < 0)
60  {
61  // Should not happen, but happens sometimes (numerical precision)
62 
63  return true;
64  }
65  double t = (-b + std::sqrt(delta)) / (2.0 * a);
66  if(t < 0.0)
67  {
68 
69  // Should not happen, but happens sometimes (numerical precision)
70  return true;
71  }
72  if(t >= 1.0)
73  {
74  // Inside the sphere
75  return false;
76  }
77 
78  if(t < 1.e-16)
79  {
80 
81  t = 0.01;
82  }
83 
84  v = v * t;
85 
86  return true;
87 }
88 
106 template< typename VertexIterator,
107  typename HalfedgeGraph,
108  typename PointMap,
109  typename MSDM2NearestMap,
110  typename MSDM2Map >
111 void
112 compute_statistics(const HalfedgeGraph &/*m_degr*/,
113  const PointMap &pm_degr,
114  const MSDM2NearestMap &nearest_pmap,
115  MSDM2Map &msdm2_map,
116  const VertexIterator p_vertex,
117  double param,
118  std::vector< double > &tab_distance1,
119  std::vector< double > &tab_distance2,
120  std::vector< CGALPoint > &tab_point1,
121  std::vector< CGALPoint > &tab_point2,
122  double radius,
123  double /*dim*/)
124 {
125  double moyenne1 = 0;
126  double variance1 = 0;
127 
128  double moyenne2 = 0;
129  double variance2 = 0;
130 
131  double covariance = 0;
132 
134 
135  double variance = radius / 2;
136  double *tab_wi1 = new double[tab_point1.size()];
137  double *tab_wi2 = new double[tab_point1.size()];
138  double somme_distance1 = 0;
139  double somme_distance2 = 0;
140  double somme_wi1 = 0;
141  double somme_wi2 = 0;
142  for(unsigned int i = 0; i < tab_point1.size(); i++) // MT
143  {
144  auto distance_pt1 = tab_point1[i] - get(pm_degr, p_vertex);
145  double dist_pt1 = sqrt(distance_pt1 * distance_pt1);
146  double wi1 = 1 / variance / sqrt(2 * 3.141592) *
147  exp(-(dist_pt1 * dist_pt1) / 2 / variance / variance);
148 
149  tab_wi1[i] = wi1;
150  somme_wi1 += wi1;
151  somme_distance1 += tab_distance1[i] * wi1;
152 
153  auto distance_pt2 = tab_point2[i] - get(nearest_pmap, p_vertex);
154  double dist_pt2 = sqrt(distance_pt2 * distance_pt2);
155 
156  double wi2 = 1 / variance / sqrt(2 * 3.141592) *
157  exp(-(dist_pt2 * dist_pt2) / 2 / variance / variance);
158  tab_wi2[i] = wi2;
159  somme_wi2 += wi2;
160  somme_distance2 += tab_distance2[i] * wi2;
161  }
162  moyenne1 = somme_distance1 / (double)somme_wi1;
163  moyenne2 = somme_distance2 / (double)somme_wi2;
164 
165 
166  for(unsigned int i = 0; i < tab_point1.size(); i++) // MT
167  {
168  variance1 += tab_wi1[i] * pow(tab_distance1[i] - moyenne1, 2);
169  variance2 += tab_wi2[i] * pow(tab_distance2[i] - moyenne2, 2);
170  covariance += tab_wi2[i] * (tab_distance1[i] - moyenne1) *
171  (tab_distance2[i] - moyenne2);
172  }
173 
174  variance1 = variance1 / somme_wi1;
175  variance2 = variance2 / somme_wi2;
176  variance1 = sqrt(variance1);
177  variance2 = sqrt(variance2);
178 
179  covariance = covariance / somme_wi1;
180 
182  // double C1=1; // MT
183  // double C2=1; // MT
184  // double C3=0.5; // MT
185  double fact1 =
186  (std::fabs(moyenne1 - moyenne2)) / (std::max(moyenne1, moyenne2) + 1);
187  double fact2 =
188  (std::fabs(variance1 - variance2)) / (std::max(variance1, variance2) + 1);
189  double fact3 = (std::fabs(variance1 * variance2 - covariance)) /
190  (variance1 * variance2 + 1);
191 
192 
193  auto tmp_msdm2 = get(msdm2_map, p_vertex);
194  tmp_msdm2 += std::pow(
195  (std::pow(fact1, 1) + std::pow(fact2, 1) + param * std::pow(fact3, 1)) /
196  (2. + param),
197  1. / 1.);
198  put(msdm2_map, p_vertex, tmp_msdm2);
199 
200 
201  delete[] tab_wi1;
202  delete[] tab_wi2;
203 }
204 
205 
216 template< typename MSDM2Map,
217  typename HalfedgeGraph,
218  typename MSDM2NearestMap,
219  typename PointMap,
220  typename TagMap,
221  typename FaceIterator =
222  typename boost::graph_traits< HalfedgeGraph >::face_iterator >
223 void
224 matching_multires_init(const HalfedgeGraph &m_poly_degrad,
225  MSDM2Map &msdm2m_degrad,
226  MSDM2NearestMap &msdm2nm_degrad,
227  TagMap &tag_map,
228  const PointMap &pm_degrad,
229  const HalfedgeGraph &m_poly_original,
230  CGAL::SM_Face_index *tab_matched_facet)
231 {
232  using AABB_Primitive =
233  typename CGAL::AABB_face_graph_triangle_primitive< HalfedgeGraph >;
234  using AABB_Traits = typename CGAL::AABB_traits< CGALKernel, AABB_Primitive >;
235  using AABB_Tree = typename CGAL::AABB_tree< AABB_Traits >;
236  //using Object_and_primitive_id = typename AABB_Tree::Object_and_primitive_id;
237  using Point_and_primitive_id = typename AABB_Tree::Point_and_primitive_id;
238 
239 
240  auto original_faces = faces(m_poly_original);
241 
242  FaceIterator orig_begin = original_faces.first;
243  FaceIterator orig_end = original_faces.second;
244 
245  AABB_Tree tree(orig_begin, orig_end, m_poly_original);
246  tree.accelerate_distance_queries();
247 
248  // Searching for the closest point and facet for each vertex
249 
250  int ind = 0;
251  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
252  {
253  put(msdm2m_degrad, vi, 0.);
254  put(tag_map, vi, ind);
255  // computes closest point and primitive id
256  Point_and_primitive_id pp =
257  tree.closest_point_and_primitive(get(pm_degrad, vi));
258  auto f_nearest = pp.second; // closest primitive id
259 
260  put(msdm2nm_degrad, vi, pp.first);
261  tab_matched_facet[ind] = f_nearest;
262 
263  ind++;
264  }
265 }
266 
267 
268 template< typename PointMap,
269  typename HalfedgeGraph,
270  typename FaceNormalMap,
271  typename VertexMSDM2Map,
272  typename GeometryTraits >
273 double
274 process_msdm2_multires(const HalfedgeGraph &m_poly_degrad,
275  const PointMap &pm_degrad,
276  const FaceNormalMap &fnm_degrad,
277  FEVV::PMapsContainer &pmaps_degrad,
278  const GeometryTraits &gt_degrad,
279  const HalfedgeGraph &m_poly_original,
280  const PointMap &pm_original,
281  const FaceNormalMap &fnm_original,
282  FEVV::PMapsContainer &/*pmaps_original*/,
283  const GeometryTraits &/*gt_original*/,
284  const int nb_level,
285  VertexMSDM2Map &msdm2_pmap,
286  const double maxdim,
287  double &msdm2_value)
288 {
289  std::vector< CGALPoint > tab_point1;
290  std::vector< CGALPoint > tab_point2;
291 
292  double somme_msdm2 = 0;
293  int nb_vertex = 0;
294 
295  auto curvature_map_degr =
296  FEVV::make_vertex_property_map< HalfedgeGraph,
298  m_poly_degrad);
299  auto curvature_map_orig =
300  FEVV::make_vertex_property_map< HalfedgeGraph,
302  m_poly_original);
303 
304  CGAL::SM_Face_index *tab_matched_facet =
305  new CGAL::SM_Face_index[size_of_vertices(m_poly_degrad)];
306 
307 
308  VertexMSDM2Map msdm2_match_pmap;
309  if(has_map(pmaps_degrad, std::string("v:msdm2_match")))
310  {
311  msdm2_match_pmap =
312  boost::any_cast< VertexMSDM2Map >(pmaps_degrad.at("v:msdm2_match"));
313  }
314  else
315  {
316  msdm2_match_pmap =
317  FEVV::make_vertex_property_map< HalfedgeGraph, double >(m_poly_degrad);
318  pmaps_degrad["v:msdm2_match"] = msdm2_match_pmap;
319  }
320 
322  vertex_msdm2_nearest_map;
323  vertex_msdm2_nearest_map msdm2_nearest_pmap;
324 
325  if(has_map(pmaps_degrad, std::string("v:msdm2_nearest")))
326  {
327  msdm2_nearest_pmap = boost::any_cast< vertex_msdm2_nearest_map >(
328  pmaps_degrad.at("v:msdm2_nearest"));
329  }
330  else
331  {
332  msdm2_nearest_pmap =
333  FEVV::make_vertex_property_map< HalfedgeGraph, CGALPoint >(
334  m_poly_degrad);
335  pmaps_degrad["v:msdm2_nearest"] = msdm2_nearest_pmap;
336  }
337 
338  typedef typename FEVV::Vertex_pmap< HalfedgeGraph, int > vertex_tag_map;
339  vertex_tag_map tag_map;
340 
341  if(has_map(pmaps_degrad, std::string("v:tag")))
342  {
343  tag_map = boost::any_cast< vertex_tag_map >(pmaps_degrad.at("v:tag"));
344  }
345  else
346  {
347  tag_map =
348  FEVV::make_vertex_property_map< HalfedgeGraph, int >(m_poly_degrad);
349  pmaps_degrad["v:tag"] = tag_map;
350  }
351 
352  matching_multires_init(m_poly_degrad,
353  msdm2_pmap,
354  msdm2_nearest_pmap,
355  tag_map,
356  pm_degrad,
357  m_poly_original,
358  tab_matched_facet);
359 
360  double radius_curvature = 0.002;
361  for(int i = 0; i < nb_level; i++)
362  {
363  double min_nrm_min_curvature_deg, max_nrm_min_curvature_deg,
364  min_nrm_max_curvature_deg, max_nrm_max_curvature_deg;
365  double min_nrm_min_curvature_org, max_nrm_min_curvature_org,
366  min_nrm_max_curvature_org, max_nrm_max_curvature_org;
367 
369  curvature_map_degr,
370  pm_degrad,
371  fnm_degrad,
372  true,
373  maxdim * radius_curvature,
374  min_nrm_min_curvature_deg,
375  max_nrm_min_curvature_deg,
376  min_nrm_max_curvature_deg,
377  max_nrm_max_curvature_deg);
378  FEVV::Filters::calculate_curvature(m_poly_original,
379  curvature_map_orig,
380  pm_original,
381  fnm_original,
382  true,
383  maxdim * radius_curvature,
384  min_nrm_min_curvature_org,
385  max_nrm_min_curvature_org,
386  min_nrm_max_curvature_org,
387  max_nrm_max_curvature_org);
388 
389  kmax_kmean(m_poly_original, curvature_map_orig, maxdim);
390  kmax_kmean(m_poly_degrad, curvature_map_degr, maxdim);
391 
392  matching_multires_update(m_poly_degrad,
393  pm_degrad,
394  m_poly_original,
395  pm_original,
396  curvature_map_orig,
397  msdm2_match_pmap,
398  tab_matched_facet);
399 
400  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
401  {
402  std::vector< double > tab_distance1;
403  std::vector< double > tab_distance2;
404 
405  tab_point1.clear();
406  tab_point2.clear();
407 
408  process_msdm2_per_vertex(m_poly_degrad,
409  pm_degrad,
410  tag_map,
411  curvature_map_degr,
412  msdm2_nearest_pmap,
413  msdm2_match_pmap,
414  gt_degrad,
415  m_poly_original,
416  pm_original,
417  vi,
418  radius_curvature * 5 * maxdim,
419  tab_distance1,
420  tab_distance2,
421  tab_point1,
422  tab_point2);
423  compute_statistics(m_poly_degrad,
424  pm_degrad,
425  msdm2_nearest_pmap,
426  msdm2_pmap,
427  vi,
428  0.5,
429  tab_distance1,
430  tab_distance2,
431  tab_point1,
432  tab_point2,
433  radius_curvature * 5 * maxdim,
434  maxdim);
435  }
436 
437 
438  radius_curvature += 0.001;
439  }
440 
441  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
442  {
443  double l_msdm2 = get(msdm2_pmap, vi) / nb_level;
444  put(msdm2_pmap, vi, l_msdm2);
445  somme_msdm2 += std::pow(l_msdm2, 3);
446  nb_vertex++;
447  }
448 
449 
450  msdm2_value = std::pow(somme_msdm2 / (double)nb_vertex, 1. / 3.);
451 
452  delete[] tab_matched_facet;
453 
454  return 0;
455 }
456 
466 template< typename PointMap,
467  typename HalfedgeGraph,
468  typename CurvatureVertexMap,
469  typename MSMD2CurvatureMatchMap >
470 void
471 matching_multires_update(const HalfedgeGraph &m_poly_degrad,
472  const PointMap &pm_degr,
473  const HalfedgeGraph &m_poly_orig,
474  const PointMap &pm_orig,
475  const CurvatureVertexMap curv_orig,
476  MSMD2CurvatureMatchMap curvmatch_pmap,
477  CGAL::SM_Face_index *tab_matched_facet)
478 {
479 
480  int ind = 0;
481  BOOST_FOREACH(Vertex_Descriptor vi, vertices(m_poly_degrad))
482  {
483  // Point3d Nearest=pVertex->match; // MT
484  auto *f_nearest = &tab_matched_facet[ind];
485 
486  ind++;
489  // we use linear interpolation using barycentric coordinates
490  auto x1 =
491  get(pm_orig, target(halfedge(*f_nearest, m_poly_orig), m_poly_orig));
492  auto x2 = get(pm_orig,
493  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
494  m_poly_orig));
495  auto x3 =
496  get(pm_orig,
497  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
498  m_poly_orig),
499  m_poly_orig));
500 
501  double l1 = sqrt((x3 - x2) * (x3 - x2));
502  double l2 = sqrt((x1 - x3) * (x1 - x3));
503  double l3 = sqrt((x1 - x2) * (x1 - x2));
504 
505  auto v1 = x1 - get(pm_degr, vi);
506  auto v2 = x2 - get(pm_degr, vi);
507  auto v3 = x3 - get(pm_degr, vi);
508 
509  double t1 = sqrt(v1 * v1);
510  double t2 = sqrt(v2 * v2);
511  double t3 = sqrt(v3 * v3);
512 
513  double p1 = (l1 + t2 + t3) / 2;
514  double p2 = (t1 + l2 + t3) / 2;
515  double p3 = (t1 + t2 + l3) / 2;
516 
517  double a1 = (p1 * (p1 - l1) * (p1 - t3) * (p1 - t2));
518  double a2 = (p2 * (p2 - l2) * (p2 - t3) * (p2 - t1));
519  double a3 = (p3 * (p3 - l3) * (p3 - t1) * (p3 - t2));
520 
521  if(a1 > 0)
522  a1 = sqrt(a1);
523  else
524  a1 = 0;
525  if(a2 > 0)
526  a2 = sqrt(a2);
527  else
528  a2 = 0;
529  if(a3 > 0)
530  a3 = sqrt(a3);
531  else
532  a3 = 0;
533 
534  double c1 =
535  get(curv_orig, target(halfedge(*f_nearest, m_poly_orig), m_poly_orig))
536  .KmaxCurv;
537  double c2 = get(curv_orig,
538  target(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
539  m_poly_orig))
540  .KmaxCurv;
541  double c3 =
542  get(curv_orig,
543  target(next(next(halfedge(*f_nearest, m_poly_orig), m_poly_orig),
544  m_poly_orig),
545  m_poly_orig))
546  .KmaxCurv;
547 
548  double curvmatch = 0.0;
549  if((a1 + a2 + a3) > 0)
550  curvmatch = (a1 * c1 + a2 * c2 + a3 * c3) / (a1 + a2 + a3);
551  else
552  curvmatch = (c1 + c2 + c3) / 3;
553  put(curvmatch_pmap, vi, curvmatch);
554  }
555 }
556 
572 template< typename VertexIterator,
573  typename HalfedgeGraph,
574  typename PointMap,
575  typename TagMap,
576  typename CurvatureVertexMap,
577  typename MSDM2NearestMap,
578  typename MSMD2CurvatureMatchMap,
579  typename GeometryTraits >
580 void
581 process_msdm2_per_vertex(const HalfedgeGraph &m_degr,
582  const PointMap &pm_degr,
583  const TagMap &tags_pmap,
584  const CurvatureVertexMap &curv_pmap,
585  const MSDM2NearestMap &nearest_pmap,
586  const MSMD2CurvatureMatchMap &curvmatch_pmap,
587  const GeometryTraits &/*gt_degr*/,
588  const HalfedgeGraph &/*m_orig*/,
589  const PointMap &/*pm_orig*/,
590  const VertexIterator &p_vertex,
591  double radius,
592  std::vector< double > &tab_distance1,
593  std::vector< double > &tab_distance2,
594  std::vector< CGALPoint > &tab_point1,
595  std::vector< CGALPoint > &tab_point2)
596 {
597 
598  std::set< int > vertices;
599 
600  std::stack< VertexIterator > s;
601 
602  auto o = get(pm_degr, p_vertex);
603 
604  s.push(p_vertex);
605  vertices.insert(get(tags_pmap, p_vertex));
606 
607 
608  tab_distance1.push_back(get(curv_pmap, p_vertex).KmaxCurv);
609  tab_point1.push_back(get(pm_degr, p_vertex));
610 
611  tab_distance2.push_back(get(curvmatch_pmap, p_vertex));
612  tab_point2.push_back(get(nearest_pmap, p_vertex));
613 
614  int nb_sommet_in_sphere = 0;
615 
616 
617  while(!s.empty())
618  {
619  VertexIterator v_it = s.top();
620  s.pop();
621  CGALPoint p = get(pm_degr, v_it);
622  auto h = halfedge(v_it, m_degr);
623  auto h0 = h;
624  do
625  {
626  auto p1 = get(pm_degr, target(h, m_degr));
627  auto p2 = get(pm_degr, target(opposite(h, m_degr), m_degr));
628 
629  auto p1m = get(nearest_pmap, target(h, m_degr));
630  auto p2m = get(nearest_pmap, target(opposite(h, m_degr), m_degr));
631 
632  auto v = (p2 - p1);
633  auto vm = (p2m - p1m);
634 
635  if(v_it == p_vertex || v * (p - o) > 0.0)
636  {
637 
638  double len_old = std::sqrt(v * v);
639  bool isect = sphere_clip_vector_msdm2(o, radius, p, v);
640  double len_edge = std::sqrt(v * v);
641 
642  nb_sommet_in_sphere++;
643 
644 
645  double weighted_curv1, weighted_curv2;
646  CGALPoint weighted_p1, weighted_p2;
647 
648  bool is_already_integrated = false;
649  if(!isect)
650  {
651 
652  auto w = target(opposite(h, m_degr), m_degr);
653  if(vertices.find(get(tags_pmap, w)) == vertices.end())
654  {
655  vertices.insert(get(tags_pmap, w));
656  s.push(w);
657  }
658  else
659  is_already_integrated = true;
660  }
661 
662  if(is_already_integrated == false)
663  {
664  if(len_old != 0 && isect)
665  {
666  weighted_curv1 =
667  (1 - len_edge / len_old) *
668  get(curv_pmap, target(h, m_degr)).KmaxCurv +
669  len_edge / len_old *
670  get(curv_pmap, target(opposite(h, m_degr), m_degr))
671  .KmaxCurv;
672  weighted_p1 = p1 + v;
673 
674  weighted_curv2 =
675  (1 - len_edge / len_old) *
676  get(curvmatch_pmap, target(h, m_degr)) +
677  len_edge / len_old *
678  get(curvmatch_pmap, target(opposite(h, m_degr), m_degr));
679  weighted_p2 = p1m + (len_edge / len_old) * vm;
680  }
681  else
682  {
683  weighted_curv1 =
684  get(curv_pmap, target(opposite(h, m_degr), m_degr)).KmaxCurv;
685  weighted_curv2 =
686  get(curvmatch_pmap, target(opposite(h, m_degr), m_degr));
687  weighted_p1 = p2;
688  weighted_p2 = p2m;
689  }
690 
691  tab_distance1.push_back(weighted_curv1);
692  tab_point1.push_back(weighted_p1);
693 
694  tab_distance2.push_back(weighted_curv2);
695  tab_point2.push_back(weighted_p2);
696  }
697  }
698  h = opposite(next(h, m_degr), m_degr);
699  } while(h != h0);
700  }
701 }
702 
710 template< typename HalfedgeGraph, typename CurvatureVertexMap >
711 void
712 kmax_kmean(const HalfedgeGraph &mesh,
713  CurvatureVertexMap &curv_vmap,
714  const double coef)
715 {
716 
717  BOOST_FOREACH(Vertex_Descriptor vi, vertices(mesh))
718  {
719  auto curv = get(curv_vmap, vi);
720  double kmax = curv.KmaxCurv * coef;
721  double kmin = std::abs(curv.KminCurv * coef);
722  curv.KmaxCurv = (kmax + kmin) / 2.;
723  put(curv_vmap, vi, curv);
724  }
725 }
726 
736 template< typename HalfedgeGraph,
737  typename GeometryTraits = FEVV::Geometry_traits< HalfedgeGraph >,
738  typename PointMap,
739  typename FaceNormalMap,
740  typename VertexMSDM2Map >
741 double
742 process_msdm2_multires(const HalfedgeGraph &m_poly_degrad,
743  const PointMap &pm_degrad,
744  const FaceNormalMap &fnm_degrad,
745  FEVV::PMapsContainer &pmaps_degrad,
746  const HalfedgeGraph &m_poly_original,
747  const PointMap &pm_original,
748  const FaceNormalMap &fnm_original,
749  FEVV::PMapsContainer &pmaps_original,
750  const int nb_level,
751  VertexMSDM2Map &msdm2_pmap,
752  double &msdm2_value)
753 {
754 
755  GeometryTraits gt_deg(m_poly_degrad);
756  GeometryTraits gt_org(m_poly_original);
757 
758  return process_msdm2_multires(m_poly_degrad,
759  pm_degrad,
760  fnm_degrad,
761  pmaps_degrad,
762  gt_deg,
763  m_poly_original,
764  pm_original,
765  fnm_original,
766  pmaps_original,
767  gt_org,
768  nb_level,
769  msdm2_pmap,
770  Filters::get_max_bb_size(m_poly_degrad,
771  pm_degrad,
772  gt_deg),
773  msdm2_value);
774 }
775 
776 } // namespace Filters
777 } // 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::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::Filters::kmax_kmean
void kmax_kmean(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: msdm2.h:712
FEVV::Filters::v_Curv
Definition: curvature.hpp:41
FEVV::Filters::compute_statistics
void compute_statistics(const HalfedgeGraph &, const PointMap &pm_degr, const MSDM2NearestMap &nearest_pmap, MSDM2Map &msdm2_map, const VertexIterator p_vertex, double param, std::vector< double > &tab_distance1, std::vector< double > &tab_distance2, std::vector< CGALPoint > &tab_point1, std::vector< CGALPoint > &tab_point2, double radius, double)
Calculates the local curvature statistics per vertex.
Definition: msdm2.h:112
FEVV::Vertex_Descriptor
MeshSurface::Vertex_index Vertex_Descriptor
Definition: cmdm.h:44
Enriched_Triangle
An enriched triangle.
Definition: boolops_polyhedra.hpp:59
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::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
FEVV::PMapsContainer
std::map< std::string, boost::any > PMapsContainer
Definition: properties.h:99
curvature.hpp
FEVV::Filters::calculate_curvature
void calculate_curvature(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.hpp:312
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
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
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::Filters::matching_multires_init
void matching_multires_init(const HalfedgeGraph &m_poly_degrad, MSDM2Map &msdm2m_degrad, MSDM2NearestMap &msdm2nm_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: msdm2.h:224
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::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::process_msdm2_multires
double process_msdm2_multires(const HalfedgeGraph &m_poly_degrad, const PointMap &pm_degrad, const FaceNormalMap &fnm_degrad, FEVV::PMapsContainer &pmaps_degrad, const GeometryTraits &gt_degrad, const HalfedgeGraph &m_poly_original, const PointMap &pm_original, const FaceNormalMap &fnm_original, FEVV::PMapsContainer &, const GeometryTraits &, const int nb_level, VertexMSDM2Map &msdm2_pmap, const double maxdim, double &msdm2_value)
Definition: msdm2.h:274
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::Filters::sphere_clip_vector_msdm2
bool sphere_clip_vector_msdm2(Point3d &o, double r, const Point3d &p, Vector &v)
Definition: msdm2.h:47
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
Point3d
EnrichedPolyhedron::Point_3 Point3d
Definition: boolops_cpolyhedron_builder.hpp:27
get_max_bb_size.hpp
FEVV::Segment
CGALKernel::Segment_3 Segment
Definition: cmdm.h:42
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
FEVV::Filters::process_msdm2_per_vertex
void process_msdm2_per_vertex(const HalfedgeGraph &m_degr, const PointMap &pm_degr, const TagMap &tags_pmap, const CurvatureVertexMap &curv_pmap, const MSDM2NearestMap &nearest_pmap, const MSMD2CurvatureMatchMap &curvmatch_pmap, const GeometryTraits &, const HalfedgeGraph &, const PointMap &, const VertexIterator &p_vertex, double radius, std::vector< double > &tab_distance1, std::vector< double > &tab_distance2, std::vector< CGALPoint > &tab_point1, std::vector< CGALPoint > &tab_point2)
Computes the local neighborhoods.
Definition: msdm2.h:581