blob: ff8bb139b7a582e3843c509b599983d175371b85 (
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
95
96
97
98
99
100
101
102
103
104
105
106
|
/* 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 (France)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program. If not, see <http://www.gnu.org/licenses/>.
*/
#ifndef EUCLIDEAN_WITNESS_COMPLEX_H_
#define EUCLIDEAN_WITNESS_COMPLEX_H_
#include <gudhi/Witness_complex.h>
#include <gudhi/Active_witness/Active_witness.h>
#include <gudhi/Kd_tree_search.h>
#include <utility>
#include <vector>
#include <list>
#include <limits>
namespace Gudhi {
namespace witness_complex {
/**
* \private
* \class Euclidean_witness_complex
* \brief Constructs (weak) 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_witness_complex
: public 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 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_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.
* @param[in] vertex Vertex handle of the point to retrieve.
*/
Point_d get_point(Vertex_handle vertex) const {
return landmarks_[vertex];
}
//@}
};
} // namespace witness_complex
} // namespace Gudhi
#endif // EUCLIDEAN_WITNESS_COMPLEX_H_
|