File: Local_eigen_analysis.h

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 (659 lines) | stat: -rw-r--r-- 24,400 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
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
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
// Copyright (c) 2017 GeometryFactory Sarl (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1/Classification/include/CGAL/Classification/Local_eigen_analysis.h $
// $Id: include/CGAL/Classification/Local_eigen_analysis.h b26b07a1242 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s)     : Simon Giraudot

#ifndef CGAL_CLASSIFICATION_LOCAL_EIGEN_ANALYSIS_H
#define CGAL_CLASSIFICATION_LOCAL_EIGEN_ANALYSIS_H

#include <CGAL/license/Classification.h>

#include <vector>
#include <memory>

#include <CGAL/Classification/compressed_float.h>
#include <CGAL/Search_traits_3.h>
#include <CGAL/Fuzzy_sphere.h>
#include <CGAL/Orthogonal_k_neighbor_search.h>
#include <CGAL/Default_diagonalize_traits.h>
#include <CGAL/centroid.h>
#include <CGAL/squared_distance_3.h>
#include <CGAL/array.h>
#include <CGAL/PCA_util.h>
#include <CGAL/boost/graph/properties.h>

#include <boost/graph/graph_traits.hpp>

#ifdef CGAL_LINKED_WITH_TBB
#include <tbb/parallel_for.h>
#include <tbb/blocked_range.h>
#include <tbb/scalable_allocator.h>
#include <mutex>

#endif // CGAL_LINKED_WITH_TBB

namespace CGAL {

namespace Classification {

  /*!
    \ingroup PkgClassificationDataStructures

    \brief Class that precomputes and stores the eigenvectors and
    eigenvalues of the covariance matrices of all points of a point
    set using a local neighborhood.

    This class can be used to compute eigen features (see
    \ref PkgClassificationFeatures) and to estimate local normal vectors
    and tangent planes.

  */
class Local_eigen_analysis
{
public:
  using Eigenvalues = std::array<float, 3>; ///< Eigenvalues (sorted in ascending order)

private:

#ifdef CGAL_LINKED_WITH_TBB
  template <typename PointRange, typename PointMap, typename NeighborQuery, typename DiagonalizeTraits>
  class Compute_eigen_values
  {
    Local_eigen_analysis& m_eigen;
    const PointRange& m_input;
    PointMap m_point_map;
    const NeighborQuery& m_neighbor_query;
    float& m_mean_range;
    std::mutex& m_mutex;

  public:

    Compute_eigen_values (Local_eigen_analysis& eigen,
                          const PointRange& input,
                          PointMap point_map,
                          const NeighborQuery& neighbor_query,
                          float& mean_range,
                          std::mutex& mutex)
      : m_eigen (eigen), m_input (input), m_point_map (point_map),
        m_neighbor_query (neighbor_query), m_mean_range (mean_range), m_mutex (mutex)
    { }

    void operator()(const tbb::blocked_range<std::size_t>& r) const
    {
      std::vector<std::size_t> neighbors;
      for (std::size_t i = r.begin(); i != r.end(); ++ i)
      {
        neighbors.clear();
        m_neighbor_query (get(m_point_map, *(m_input.begin()+i)), std::back_inserter (neighbors));

        std::vector<typename PointMap::value_type> neighbor_points;
        neighbor_points.reserve(neighbors.size());
        for (std::size_t j = 0; j < neighbors.size(); ++ j)
          neighbor_points.push_back (get(m_point_map, *(m_input.begin()+neighbors[j])));

        m_mutex.lock();
        m_mean_range += float(CGAL::sqrt
                              (CGAL::squared_distance (get(m_point_map, *(m_input.begin() + i)),
                                                       get(m_point_map, *(m_input.begin() + neighbors.back())))));
        m_mutex.unlock();

        m_eigen.compute<typename PointMap::value_type,
                        DiagonalizeTraits> (i, get(m_point_map, *(m_input.begin()+i)), neighbor_points);
      }
    }

  };

  template <typename FaceListGraph, typename NeighborQuery, typename DiagonalizeTraits>
  class Compute_eigen_values_graph
  {
    using face_descriptor = typename boost::graph_traits<FaceListGraph>::face_descriptor;
    using face_index = typename boost::property_map<FaceListGraph, CGAL::face_index_t>::type::value_type;
    using face_iterator = typename boost::graph_traits<FaceListGraph>::face_iterator;

    Local_eigen_analysis& m_eigen;
    const FaceListGraph& m_input;
    const NeighborQuery& m_neighbor_query;
    float& m_mean_range;
    std::mutex& m_mutex;

  public:

    Compute_eigen_values_graph (Local_eigen_analysis& eigen,
                                const FaceListGraph& input,
                                const NeighborQuery& neighbor_query,
                                float& mean_range,
                                std::mutex& mutex)
      : m_eigen (eigen), m_input (input),
        m_neighbor_query (neighbor_query), m_mean_range (mean_range), m_mutex (mutex)
    { }

    void operator()(const tbb::blocked_range<std::size_t>& r) const
    {
      face_iterator begin = faces(m_input).first;
      for (std::size_t i = r.begin(); i != r.end(); ++ i)
      {
        face_descriptor fd = *(begin + i);
        std::vector<face_index> neighbors;
        m_neighbor_query (fd, std::back_inserter (neighbors));

        m_mutex.lock();
        m_mean_range += m_eigen.face_radius(fd, m_input);
        m_mutex.unlock();

        m_eigen.compute_triangles<FaceListGraph, DiagonalizeTraits>
          (m_input, fd, neighbors);
      }
    }

  };
#endif

  template <typename ClusterRange, typename DiagonalizeTraits>
  class Compute_clusters_eigen_values
  {
    Local_eigen_analysis& m_eigen;
    const ClusterRange& m_input;

  public:

    Compute_clusters_eigen_values (Local_eigen_analysis& eigen,
                                   const ClusterRange& input)
      : m_eigen (eigen), m_input (input)
    { }

#ifdef CGAL_LINKED_WITH_TBB
    void operator()(const tbb::blocked_range<std::size_t>& r) const
    {
      for (std::size_t i = r.begin(); i != r.end(); ++ i)
        apply (i);
    }
#endif

    inline void apply (std::size_t i) const
    {
      using Cluster = typename ClusterRange::value_type;
      using Item = typename Cluster::Item;
      const Cluster& cluster = m_input[i];

      std::vector<typename ClusterRange::value_type::Item> points;
      for (std::size_t j = 0; j < cluster.size(); ++ j)
        points.push_back (cluster[j]);

      m_eigen.compute<Item, DiagonalizeTraits> (i, Item(0.,0.,0.), points);
    }

  };

  using float3 = std::array<float, 3>;
  using float2 = std::array<float, 2>;
  using cfloat2 = std::array<compressed_float, 2>;

  struct Content
  {
    std::vector<cfloat2> eigenvalues;
    std::vector<float3> centroids;
    std::vector<float2> smallest_eigenvectors;
    float mean_range;
  };

  std::shared_ptr<Content> m_content; // To avoid copies with named constructors

public:

  /// \cond SKIP_IN_MANUAL
  Local_eigen_analysis () { }
  /// \endcond

  /// \name Named Constructors
  /// @{

  /*!
    \brief computes the local eigen analysis of an input point set
    based on a local neighborhood.

    \tparam PointRange model of `ConstRange`. Its iterator type is
    `RandomAccessIterator` and its value type is the key type of
    `PointMap`.
    \tparam PointMap model of `ReadablePropertyMap` whose key
    type is the value type of the iterator of `PointRange` and value type
    is `CGAL::Point_3`.
    \tparam NeighborQuery model of `NeighborQuery`
    \tparam ConcurrencyTag enables sequential versus parallel
    algorithm. Possible values are `Parallel_tag` (default value if \cgal
    is linked with TBB) or `Sequential_tag` (default value otherwise).
    \tparam DiagonalizeTraits model of `DiagonalizeTraits` used for
    matrix diagonalization. It can be omitted if Eigen 3 (or greater)
    is available and `CGAL_EIGEN3_ENABLED` is defined. In that case,
    an overload using `Eigen_diagonalize_traits` is provided.

    \param input point range.
    \param point_map property map to access the input points.
    \param neighbor_query object used to access neighborhoods of points.
  */
  template <typename PointRange,
            typename PointMap,
            typename NeighborQuery,
#if defined(DOXYGEN_RUNNING)
            typename ConcurrencyTag,
#else
            typename ConcurrencyTag = CGAL::Parallel_if_available_tag,
#endif
#if defined(DOXYGEN_RUNNING)
            typename DiagonalizeTraits>
#else
            typename DiagonalizeTraits = CGAL::Default_diagonalize_traits<float, 3> >
#endif
  static Local_eigen_analysis create_from_point_set(const PointRange& input,
                                                    PointMap point_map,
                                                    const NeighborQuery& neighbor_query,
                                                    const ConcurrencyTag& = ConcurrencyTag(),
                                                    const DiagonalizeTraits& = DiagonalizeTraits())
  {
    Local_eigen_analysis out;
    out.m_content = std::make_shared<Content>();
    out.m_content->eigenvalues.resize (input.size());
    out.m_content->centroids.resize (input.size());
    out.m_content->smallest_eigenvectors.resize (input.size());

    out.m_content->mean_range = 0.;

#ifndef CGAL_LINKED_WITH_TBB
    static_assert (!(std::is_convertible<ConcurrencyTag, Parallel_tag>::value),
                               "Parallel_tag is enabled but TBB is unavailable.");
#else
    if (std::is_convertible<ConcurrencyTag,Parallel_tag>::value)
    {
      std::mutex mutex;
      Compute_eigen_values<PointRange, PointMap, NeighborQuery, DiagonalizeTraits>
        f(out, input, point_map, neighbor_query, out.m_content->mean_range, mutex);
      tbb::parallel_for(tbb::blocked_range<size_t>(0, input.size ()), f);
    }
    else
#endif
    {
      for (std::size_t i = 0; i < input.size(); i++)
      {
        std::vector<std::size_t> neighbors;
        neighbor_query (get(point_map, *(input.begin()+i)), std::back_inserter (neighbors));

        std::vector<typename PointMap::value_type> neighbor_points;
        for (std::size_t j = 0; j < neighbors.size(); ++ j)
          neighbor_points.push_back (get(point_map, *(input.begin()+neighbors[j])));

        out.m_content->mean_range += float(CGAL::sqrt (CGAL::squared_distance
                                          (get(point_map, *(input.begin() + i)),
                                           get(point_map, *(input.begin() + neighbors.back())))));

        out.compute<typename PointMap::value_type, DiagonalizeTraits>
          (i, get(point_map, *(input.begin()+i)), neighbor_points);
      }
    }
    out.m_content->mean_range /= input.size();

    return out;
  }


  /*!
    \brief computes the local eigen analysis of an input face graph
    based on a local neighborhood.

    \tparam FaceListGraph model of `FaceListGraph`.
    \tparam NeighborQuery model of `NeighborQuery`
    \tparam ConcurrencyTag enables sequential versus parallel
    algorithm. Possible values are `Parallel_tag` (default value if \cgal
    is linked with TBB) or `Sequential_tag` (default value otherwise).
    \tparam DiagonalizeTraits model of `DiagonalizeTraits` used for
    matrix diagonalization. It can be omitted: if Eigen 3 (or greater)
    is available and `CGAL_EIGEN3_ENABLED` is defined then an overload
    using `Eigen_diagonalize_traits` is provided. Otherwise, the
    internal implementation `Diagonalize_traits` is used.

    \param input face graph.
    \param neighbor_query object used to access neighborhoods of points.
  */
  template <typename FaceListGraph,
            typename NeighborQuery,
#if defined(DOXYGEN_RUNNING)
            typename ConcurrencyTag,
#else
            typename ConcurrencyTag = CGAL::Parallel_if_available_tag,
#endif
#if defined(DOXYGEN_RUNNING)
            typename DiagonalizeTraits>
#else
            typename DiagonalizeTraits = CGAL::Default_diagonalize_traits<float, 3> >
#endif
  static Local_eigen_analysis create_from_face_graph (const FaceListGraph& input,
                                                      const NeighborQuery& neighbor_query,
                                                      const ConcurrencyTag& = ConcurrencyTag(),
                                                      const DiagonalizeTraits& = DiagonalizeTraits())
  {
    using face_descriptor = typename boost::graph_traits<FaceListGraph>::face_descriptor;
    using face_iterator = typename boost::graph_traits<FaceListGraph>::face_iterator;
    using Face_range = typename CGAL::Iterator_range<face_iterator>;
    using face_index = typename boost::property_map<FaceListGraph, CGAL::face_index_t>::type::value_type;

    Local_eigen_analysis out;
    out.m_content = std::make_shared<Content>();

    Face_range range (faces(input));

    out.m_content->eigenvalues.resize (range.size());
    out.m_content->centroids.resize (range.size());
    out.m_content->smallest_eigenvectors.resize (range.size());

    out.m_content->mean_range = 0.;

#ifndef CGAL_LINKED_WITH_TBB
    static_assert (!(std::is_convertible<ConcurrencyTag, Parallel_tag>::value),
                               "Parallel_tag is enabled but TBB is unavailable.");
#else
    if (std::is_convertible<ConcurrencyTag,Parallel_tag>::value)
    {
      std::mutex mutex;
      Compute_eigen_values_graph<FaceListGraph, NeighborQuery, DiagonalizeTraits>
        f(out, input, neighbor_query, out.m_content->mean_range, mutex);

      tbb::parallel_for(tbb::blocked_range<std::size_t>(0, range.size()), f);
    }
    else
#endif
    {
      for(face_descriptor fd : range)
      {
        std::vector<face_index> neighbors;
        neighbor_query (fd, std::back_inserter (neighbors));

        out.m_content->mean_range += out.face_radius(fd, input);

        out.compute_triangles<FaceListGraph, DiagonalizeTraits>
          (input, fd, neighbors);

      }
    }

    out.m_content->mean_range /= range.size();
    return out;
  }

  /*!
    \brief computes the local eigen analysis of an input set of point
    clusters based on a local neighborhood.

    \tparam ClusterRange model of `ConstRange`. Its iterator type is
    `RandomAccessIterator` and its value type is the key type of
    `PointMap`.
    \tparam ConcurrencyTag enables sequential versus parallel
    algorithm. Possible values are `Parallel_tag` (default value if \cgal
    is linked with TBB) or `Sequential_tag` (default value otherwise).
    \tparam DiagonalizeTraits model of `DiagonalizeTraits` used for
    matrix diagonalization. It can be omitted: if Eigen 3 (or greater)
    is available and `CGAL_EIGEN3_ENABLED` is defined then an overload
    using `Eigen_diagonalize_traits` is provided. Otherwise, the
    internal implementation `Diagonalize_traits` is used.

    \param input cluster range.
  */
  template <typename ClusterRange,
#if defined(DOXYGEN_RUNNING)
            typename ConcurrencyTag,
#else
            typename ConcurrencyTag = CGAL::Parallel_if_available_tag,
#endif
#if defined(DOXYGEN_RUNNING)
            typename DiagonalizeTraits>
#else
            typename DiagonalizeTraits = CGAL::Default_diagonalize_traits<float, 3> >
#endif
  static Local_eigen_analysis create_from_point_clusters (const ClusterRange& input,
                                                          const ConcurrencyTag& = ConcurrencyTag(),
                                                          const DiagonalizeTraits& = DiagonalizeTraits())
  {
    Local_eigen_analysis out;
    out.m_content = std::make_shared<Content>();

    out.m_content->eigenvalues.resize (input.size());
    out.m_content->centroids.resize (input.size());
    out.m_content->smallest_eigenvectors.resize (input.size());

    out.m_content->mean_range = 0.;

    Compute_clusters_eigen_values<ClusterRange, DiagonalizeTraits>
    f(out, input);


#ifndef CGAL_LINKED_WITH_TBB
    static_assert (!(std::is_convertible<ConcurrencyTag, Parallel_tag>::value),
                               "Parallel_tag is enabled but TBB is unavailable.");
#else
    if (std::is_convertible<ConcurrencyTag,Parallel_tag>::value)
    {
      tbb::parallel_for(tbb::blocked_range<size_t>(0, input.size ()), f);
    }
    else
#endif
    {
      for (std::size_t i = 0; i < input.size(); ++ i)
        f.apply (i);
    }
    return out;
  }


  /// @}

  /// \name Analysis
  /// @{

  /*!
    \brief returns the estimated unoriented normal vector of the point at position `index`.
    \tparam GeomTraits model of \cgal Kernel.
  */
  template <typename GeomTraits>
  typename GeomTraits::Vector_3 normal_vector (std::size_t index) const
  {
    return typename GeomTraits::Vector_3(double(m_content->smallest_eigenvectors[index][0]),
                                         double(m_content->smallest_eigenvectors[index][1]),
                                         double(1. - (m_content->smallest_eigenvectors[index][0] +
                                                      m_content->smallest_eigenvectors[index][1])));
  }

  /*!
    \brief returns the estimated local tangent plane of the point at position `index`.
    \tparam GeomTraits model of \cgal Kernel.
  */
  template <typename GeomTraits>
  typename GeomTraits::Plane_3 plane (std::size_t index) const
  {
    return typename GeomTraits::Plane_3
      (typename GeomTraits::Point_3 (double(m_content->centroids[index][0]),
                                      double(m_content->centroids[index][1]),
                                      double(m_content->centroids[index][2])),
       normal_vector<GeomTraits>(index));
  }

  /*!
    \brief returns the normalized eigenvalues of the point at position `index`.
  */
  Eigenvalues eigenvalue (std::size_t index) const
  {
    const cfloat2& uc = m_content->eigenvalues[index];
    Eigenvalues out;
    out[1] = decompress_float(uc[0]);
    out[2] = decompress_float(uc[1]);
    out[0] = 1.f - (out[1] + out[2]);
    return out;
  }

  /// @}

  /// \cond SKIP_IN_MANUAL
  float mean_range() const { return m_content->mean_range; }
  /// \endcond

private:

  template <typename FaceListGraph>
  float face_radius (typename boost::graph_traits<FaceListGraph>::face_descriptor& fd,
                     const FaceListGraph& g)
  {
    using halfedge_descriptor = typename boost::graph_traits<FaceListGraph>::halfedge_descriptor;

    float out = 0.f;
    for(halfedge_descriptor hd : halfedges_around_face(halfedge(fd, g), g))
    {
      out = (std::max)(out,
                       float(CGAL::squared_distance (get(get (CGAL::vertex_point, g), source(hd,g)),
                                                     get(get (CGAL::vertex_point, g), target(hd,g)))));
    }
    return out;
  }

  template <typename Point, typename DiagonalizeTraits>
  void compute (std::size_t index, const Point& query, std::vector<Point>& neighbor_points)
  {
    using Vector = typename Kernel_traits<Point>::Kernel::Vector_3;

    if (neighbor_points.size() == 0)
    {
      m_content->eigenvalues[index] = make_array (compressed_float(0), compressed_float(0));
      m_content->centroids[index] = make_array(float(query.x()), float(query.y()), float(query.z()) );
      m_content->smallest_eigenvectors[index] = make_array( 0.f, 0.f );
      return;
    }

    Point centroid = CGAL::centroid (neighbor_points.begin(), neighbor_points.end());
    m_content->centroids[index] = make_array( float(centroid.x()), float(centroid.y()), float(centroid.z()) );

    std::array<float, 6> covariance = make_array( 0.f, 0.f, 0.f, 0.f, 0.f, 0.f );

    for (std::size_t i = 0; i < neighbor_points.size(); ++ i)
    {
      Vector d = neighbor_points[i] - centroid;
      covariance[0] += float(d.x () * d.x ());
      covariance[1] += float(d.x () * d.y ());
      covariance[2] += float(d.x () * d.z ());
      covariance[3] += float(d.y () * d.y ());
      covariance[4] += float(d.y () * d.z ());
      covariance[5] += float(d.z () * d.z ());
    }

    std::array<float, 3> evalues = make_array( 0.f, 0.f, 0.f );
    std::array<float, 9> evectors = make_array( 0.f, 0.f, 0.f,
                                               0.f, 0.f, 0.f,
                                               0.f, 0.f, 0.f );

    DiagonalizeTraits::diagonalize_selfadjoint_covariance_matrix
      (covariance, evalues, evectors);

    // Normalize
    float sum = evalues[0] + evalues[1] + evalues[2];
    if (sum > 0.f)
      for (std::size_t i = 0; i < 3; ++ i)
        evalues[i] = evalues[i] / sum;

    m_content->eigenvalues[index] = make_array(compress_float (evalues[1]),
                                               compress_float (evalues[2]));

    sum = evectors[0] + evectors[1] + evectors[2];
    if (sum > 0.f)
      for (std::size_t i = 0; i < 3; ++ i)
        evectors[i] = evectors[i] / sum;
    m_content->smallest_eigenvectors[index] = make_array( float(evectors[0]), float(evectors[1]) );
  }

  template <typename FaceListGraph, typename DiagonalizeTraits>
  void compute_triangles (const FaceListGraph& g,
                          typename boost::graph_traits<FaceListGraph>::face_descriptor& query,
                          std::vector<typename boost::property_map<FaceListGraph, CGAL::face_index_t>::type::value_type>& neighbor_faces)
  {
    using Point = typename boost::property_map<FaceListGraph, boost::vertex_point_t>::type::value_type;
    using Kernel = typename Kernel_traits<Point>::Kernel;
    using Triangle = typename Kernel::Triangle_3;

    using face_descriptor = typename boost::graph_traits<FaceListGraph>::face_descriptor;
    using face_iterator = typename boost::graph_traits<FaceListGraph>::face_iterator;

    if (neighbor_faces.size() == 0)
    {
      m_content->eigenvalues[get(get(CGAL::face_index,g), query)]
        = make_array(compressed_float(0), compressed_float(0));

      std::array<Triangle,1> tr
        = {{ Triangle (get(get (CGAL::vertex_point, g), target(halfedge(query, g), g)),
                       get(get (CGAL::vertex_point, g), target(next(halfedge(query, g), g), g)),
                       get(get (CGAL::vertex_point, g), target(next(next(halfedge(query, g), g), g), g))) }};
      Point c = CGAL::centroid(tr.begin(),
                               tr.end(), Kernel(), CGAL::Dimension_tag<2>());

      m_content->centroids[get(get(CGAL::face_index,g), query)] = {{ float(c.x()), float(c.y()), float(c.z()) }};

      m_content->smallest_eigenvectors[get(get(CGAL::face_index,g), query)] = {{ 0.f, 0.f }};
      return;
    }

    std::vector<Triangle> triangles;
    triangles.reserve(neighbor_faces.size());

    face_iterator begin = faces(g).first;
    for (std::size_t i = 0; i < neighbor_faces.size(); ++ i)
    {
      face_descriptor fd = *(begin + std::size_t(neighbor_faces[i]));
      triangles.push_back
        (Triangle (get(get (CGAL::vertex_point, g), target(halfedge(fd, g), g)),
                   get(get (CGAL::vertex_point, g), target(next(halfedge(fd, g), g), g)),
                   get(get (CGAL::vertex_point, g), target(next(next(halfedge(fd, g), g), g), g))));
    }

    std::array<float, 6> covariance = {{ 0.f, 0.f, 0.f, 0.f, 0.f, 0.f }};
    Point c = CGAL::centroid(triangles.begin(),
                             triangles.end(), Kernel(), CGAL::Dimension_tag<2>());

    CGAL::internal::assemble_covariance_matrix_3 (triangles.begin(), triangles.end(), covariance,
                                                  c, Kernel(), (Triangle*)nullptr, CGAL::Dimension_tag<2>(),
                                                  DiagonalizeTraits());

    m_content->centroids[get(get(CGAL::face_index,g), query)] = {{ float(c.x()), float(c.y()), float(c.z()) }};

    std::array<float, 3> evalues = {{ 0.f, 0.f, 0.f }};
    std::array<float, 9> evectors = {{ 0.f, 0.f, 0.f,
                                               0.f, 0.f, 0.f,
                                               0.f, 0.f, 0.f }};

    DiagonalizeTraits::diagonalize_selfadjoint_covariance_matrix
      (covariance, evalues, evectors);

    // Normalize
    float sum = evalues[0] + evalues[1] + evalues[2];
    if (sum > 0.f)
      for (std::size_t i = 0; i < 3; ++ i)
        evalues[i] = evalues[i] / sum;

    m_content->eigenvalues[get(get(CGAL::face_index,g), query)]
      = make_array(compress_float (evalues[1]),
                   compress_float (evalues[2]));

    sum = evectors[0] + evectors[1] + evectors[2];
    if (sum > 0.f)
      for (std::size_t i = 0; i < 3; ++ i)
        evectors[i] = evectors[i] / sum;
    m_content->smallest_eigenvectors[get(get(CGAL::face_index,g), query)] = {{ float(evectors[0]), float(evectors[1]), }};
  }

};


}

}


#endif // CGAL_CLASSIFICATION_LOCAL_EIGEN_ANALYSIS_H