summaryrefslogtreecommitdiff
path: root/include/gudhi/Tangential_complex/utilities.h
blob: 2dd4611896829a4cd93c1617da4ec5604d8d9c3f (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
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
/*    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):       Clement Jamin
 *
 *    Copyright (C) 2016 Inria
 *
 *    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 TANGENTIAL_COMPLEX_UTILITIES_H_
#define TANGENTIAL_COMPLEX_UTILITIES_H_

#include <CGAL/Dimension.h>
#include <CGAL/Combination_enumerator.h>
#include <CGAL/IO/Triangulation_off_ostream.h>

#include <boost/container/flat_set.hpp>

#include <Eigen/Core>
#include <Eigen/Eigen>

#include <set>
#include <vector>
#include <array>
#include <fstream>
#include <atomic>
#include <cmath>  // for std::sqrt

namespace Gudhi {
namespace tangential_complex {
namespace internal {

// Provides copy constructors to std::atomic so that
// it can be used in a vector
template <typename T>
struct Atomic_wrapper
: public std::atomic<T> {
  typedef std::atomic<T> Base;

  Atomic_wrapper() { }

  Atomic_wrapper(const T &t) : Base(t) { }

  Atomic_wrapper(const std::atomic<T> &a) : Base(a.load()) { }

  Atomic_wrapper(const Atomic_wrapper &other) : Base(other.load()) { }

  Atomic_wrapper &operator=(const T &other) {
    Base::store(other);
    return *this;
  }

  Atomic_wrapper &operator=(const std::atomic<T> &other) {
    Base::store(other.load());
    return *this;
  }

  Atomic_wrapper &operator=(const Atomic_wrapper &other) {
    Base::store(other.load());
    return *this;
  }
};

// Modifies v in-place
template <typename K>
typename K::Vector_d& normalize_vector(typename K::Vector_d& v,
                                       K const& k) {
  v = k.scaled_vector_d_object()(
                                 v, typename K::FT(1) / std::sqrt(k.squared_length_d_object()(v)));
  return v;
}

template<typename Kernel>
struct Basis {
  typedef typename Kernel::FT FT;
  typedef typename Kernel::Point_d Point;
  typedef typename Kernel::Vector_d Vector;
  typedef typename std::vector<Vector>::const_iterator const_iterator;

  std::size_t m_origin;
  std::vector<Vector> m_vectors;

  std::size_t origin() const {
    return m_origin;
  }

  void set_origin(std::size_t o) {
    m_origin = o;
  }

  const_iterator begin() const {
    return m_vectors.begin();
  }

  const_iterator end() const {
    return m_vectors.end();
  }

  std::size_t size() const {
    return m_vectors.size();
  }

  Vector& operator[](const std::size_t i) {
    return m_vectors[i];
  }

  const Vector& operator[](const std::size_t i) const {
    return m_vectors[i];
  }

  void push_back(const Vector& v) {
    m_vectors.push_back(v);
  }

  void reserve(const std::size_t s) {
    m_vectors.reserve(s);
  }

  Basis() { }

  Basis(std::size_t origin) : m_origin(origin) { }

  Basis(std::size_t origin, const std::vector<Vector>& vectors)
      : m_origin(origin), m_vectors(vectors) { }

  int dimension() const {
    return static_cast<int> (m_vectors.size());
  }
};

// 1st line: number of points
// Then one point per line
template <typename Kernel, typename Point_range>
std::ostream &export_point_set(
                               Kernel const& k,
                               Point_range const& points,
                               std::ostream & os,
                               const char *coord_separator = " ") {
  // Kernel functors
  typename Kernel::Construct_cartesian_const_iterator_d ccci =
      k.construct_cartesian_const_iterator_d_object();

  os << points.size() << "\n";

  typename Point_range::const_iterator it_p = points.begin();
  typename Point_range::const_iterator it_p_end = points.end();
  // For each point p
  for (; it_p != it_p_end; ++it_p) {
    for (auto it = ccci(*it_p); it != ccci(*it_p, 0); ++it)
      os << CGAL::to_double(*it) << coord_separator;

    os << "\n";
  }

  return os;
}

// Compute all the k-combinations of elements
// Output_iterator::value_type must be
// boost::container::flat_set<std::size_t>
template <typename Elements_container, typename Output_iterator>
void combinations(const Elements_container elements, int k,
                  Output_iterator combinations) {
  std::size_t n = elements.size();
  std::vector<bool> booleans(n, false);
  std::fill(booleans.begin() + n - k, booleans.end(), true);
  do {
    boost::container::flat_set<std::size_t> combination;
    typename Elements_container::const_iterator it_elt = elements.begin();
    for (std::size_t i = 0; i < n; ++i, ++it_elt) {
      if (booleans[i])
        combination.insert(*it_elt);
    }
    *combinations++ = combination;
  } while (std::next_permutation(booleans.begin(), booleans.end()));
}

}  // namespace internal
}  // namespace tangential_complex
}  // namespace Gudhi

#endif  // TANGENTIAL_COMPLEX_UTILITIES_H_