summaryrefslogtreecommitdiff
path: root/src/Witness_complex/include/gudhi/Euclidean_strong_witness_complex.h
blob: d2bf00cefda3f87d5d536dcb6776232158690a37 (plain)
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
/*    This file is part of the Gudhi Library. The Gudhi library
 *    (Geometric Understanding in Higher Dimensions) is a generic C++
 *    library for computational topology.
 *
 *    Author(s):       Siargey Kachanovich
 *
 *    Copyright (C) 2015 Inria
 *
 *    Modification(s):
 *      - YYYY/MM Author: Description of the modification
 */

#ifndef EUCLIDEAN_STRONG_WITNESS_COMPLEX_H_
#define EUCLIDEAN_STRONG_WITNESS_COMPLEX_H_

#include <gudhi/Strong_witness_complex.h>
#include <gudhi/Active_witness/Active_witness.h>
#include <gudhi/Kd_tree_search.h>

#include <utility>
#include <vector>

namespace Gudhi {

namespace witness_complex {

/**
 *  \private
 * \class Euclidean_strong_witness_complex
 * \brief Constructs strong witness complex for given sets of witnesses and landmarks in Euclidean space.
 * \ingroup witness_complex
 *
 * \tparam Kernel_ requires a <a target="_blank"
 * href="http://doc.cgal.org/latest/Kernel_d/classCGAL_1_1Epick__d.html">CGAL::Epick_d</a> class.
 */
template< class Kernel_ >
class Euclidean_strong_witness_complex
    : public Strong_witness_complex<std::vector<typename Gudhi::spatial_searching::Kd_tree_search<Kernel_,
                                                                                                  std::vector<typename Kernel_::Point_d>>::INS_range>> {
 private:
  typedef Kernel_                                                                      K;
  typedef typename K::Point_d                                                          Point_d;
  typedef std::vector<Point_d>                                                         Point_range;
  typedef Gudhi::spatial_searching::Kd_tree_search<Kernel_, Point_range>               Kd_tree;
  typedef typename Kd_tree::INS_range                                                  Nearest_landmark_range;
  typedef typename std::vector<Nearest_landmark_range>                                 Nearest_landmark_table;

  typedef typename Nearest_landmark_range::Point_with_transformed_distance Id_distance_pair;
  typedef typename Id_distance_pair::first_type Landmark_id;
  typedef Landmark_id Vertex_handle;

 private:
  Point_range                         landmarks_;
  Kd_tree                             landmark_tree_;
  using Strong_witness_complex<Nearest_landmark_table>::nearest_landmark_table_;

 public:
  /////////////////////////////////////////////////////////////////////////////////////////////////////////////////////
  /* @name Constructor
   */

  //@{

  /**
   *  \brief Initializes member variables before constructing simplicial complex.
   *  \details Records landmarks from the range 'landmarks' into a 
   *           table internally, as well as witnesses from the range 'witnesses'.
   *           Both ranges should have value_type Kernel_::Point_d.
   */
  template< typename LandmarkRange,
            typename WitnessRange >
  Euclidean_strong_witness_complex(const LandmarkRange & landmarks,
                                   const WitnessRange &  witnesses)
    : landmarks_(std::begin(landmarks), std::end(landmarks)), landmark_tree_(landmarks_) {
    nearest_landmark_table_.reserve(boost::size(witnesses));
    for (auto w : witnesses)
      nearest_landmark_table_.push_back(landmark_tree_.incremental_nearest_neighbors(w));
  }

  /** \brief Returns the point corresponding to the given vertex.
   */
  template <typename Vertex_handle>
  Point_d get_point(Vertex_handle vertex) const {
    return landmarks_[vertex];
  }

  //@}
};

}  // namespace witness_complex

}  // namespace Gudhi

#endif  // EUCLIDEAN_STRONG_WITNESS_COMPLEX_H_