summaryrefslogtreecommitdiff
path: root/src/Alpha_complex/test/Alpha_complex_unit_test.cpp
blob: 9530314c8e4aaf004490ca6ccb17349a99f4c18a (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
/*    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):       Vincent Rouvreau
 *
 *    Copyright (C) 2015  INRIA Saclay (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/>.
 */

#define BOOST_TEST_MODULE alpha_complex
#include <boost/test/included/unit_test.hpp>
#include <boost/system/error_code.hpp>
#include <boost/chrono/thread_clock.hpp>
// to construct a Delaunay_triangulation from a OFF file
#include "gudhi/Delaunay_triangulation_off_io.h"
#include "gudhi/Alpha_complex.h"

#include <CGAL/Delaunay_triangulation.h>
#include <CGAL/Epick_d.h>

#include <cmath> // float comparison
#include <limits>

// Use dynamic_dimension_tag for the user to be able to set dimension
typedef CGAL::Epick_d< CGAL::Dynamic_dimension_tag > Kernel;
typedef Kernel::Point_d Point;
typedef std::vector<Point> Vector_of_points;
// The triangulation uses the default instantiation of the TriangulationDataStructure template parameter

BOOST_AUTO_TEST_CASE(S4_100_OFF_file) {
  // ----------------------------------------------------------------------------
  //
  // Init of an alpha-complex from a OFF file
  //
  // ----------------------------------------------------------------------------
  std::string off_file_name("S4_100.off");
  std::cout << "========== OFF FILE NAME = " << off_file_name << " ==========" << std::endl;

  Gudhi::alphacomplex::Alpha_complex alpha_complex_from_file(off_file_name);

  const int DIMENSION = 4;
  std::cout << "alpha_complex_from_file.dimension()=" << alpha_complex_from_file.dimension() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.dimension() == DIMENSION);

  const int NUMBER_OF_VERTICES = 100;
  std::cout << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.num_vertices() == NUMBER_OF_VERTICES);

  const int NUMBER_OF_SIMPLICES = 6879;
  std::cout << "alpha_complex_from_file.num_simplices()=" << alpha_complex_from_file.num_simplices() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.num_simplices() == NUMBER_OF_SIMPLICES);

}

BOOST_AUTO_TEST_CASE(S8_10_OFF_file) {
  // ----------------------------------------------------------------------------
  //
  // Init of an alpha-complex from a OFF file
  //
  // ----------------------------------------------------------------------------
  std::string off_file_name("S8_10.off");
  std::cout << "========== OFF FILE NAME = " << off_file_name << " ==========" << std::endl;

  Gudhi::alphacomplex::Alpha_complex alpha_complex_from_file(off_file_name);

  const int DIMENSION = 8;
  std::cout << "alpha_complex_from_file.dimension()=" << alpha_complex_from_file.dimension() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.dimension() == DIMENSION);

  const int NUMBER_OF_VERTICES = 10;
  std::cout << "alpha_complex_from_file.num_vertices()=" << alpha_complex_from_file.num_vertices() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.num_vertices() == NUMBER_OF_VERTICES);

  const int NUMBER_OF_SIMPLICES = 1007;
  std::cout << "alpha_complex_from_file.num_simplices()=" << alpha_complex_from_file.num_simplices() << std::endl;
  BOOST_CHECK(alpha_complex_from_file.num_simplices() == NUMBER_OF_SIMPLICES);
}

bool are_almost_the_same(float a, float b) {
  return std::fabs(a - b) < std::numeric_limits<float>::epsilon();
}

BOOST_AUTO_TEST_CASE(Alpha_complex_from_points) {

  // ----------------------------------------------------------------------------
  // Init of a list of points
  // ----------------------------------------------------------------------------
  Vector_of_points points;
  std::vector<double> coords;

  coords.clear();
  coords.push_back(0.0);
  coords.push_back(0.0);
  coords.push_back(0.0);
  coords.push_back(1.0);
  points.push_back(Point(coords.begin(), coords.end()));
  coords.clear();
  coords.push_back(0.0);
  coords.push_back(0.0);
  coords.push_back(1.0);
  coords.push_back(0.0);
  points.push_back(Point(coords.begin(), coords.end()));
  coords.clear();
  coords.push_back(0.0);
  coords.push_back(1.0);
  coords.push_back(0.0);
  coords.push_back(0.0);
  points.push_back(Point(coords.begin(), coords.end()));
  coords.clear();
  coords.push_back(1.0);
  coords.push_back(0.0);
  coords.push_back(0.0);
  coords.push_back(0.0);
  points.push_back(Point(coords.begin(), coords.end()));

  // ----------------------------------------------------------------------------
  // Init of an alpha complex from the list of points
  // ----------------------------------------------------------------------------
  Gudhi::alphacomplex::Alpha_complex alpha_complex_from_points(3, points.size(), points.begin(), points.end());

  std::cout << "========== Alpha_complex_from_points ==========" << std::endl;

  std::cout << "alpha_complex_from_points.dimension()=" << alpha_complex_from_points.dimension() << std::endl;
  BOOST_CHECK(alpha_complex_from_points.dimension() == 3);
  std::cout << "alpha_complex_from_points.num_simplices()=" << alpha_complex_from_points.num_simplices() << std::endl;
  BOOST_CHECK(alpha_complex_from_points.num_simplices() == 15);
  std::cout << "alpha_complex_from_points.num_vertices()=" << alpha_complex_from_points.num_vertices() << std::endl;
  BOOST_CHECK(alpha_complex_from_points.num_vertices() == 4);

  for (auto f_simplex : alpha_complex_from_points.filtration_simplex_range()) {
    switch (alpha_complex_from_points.dimension(f_simplex)) {
      case 0:
        BOOST_CHECK(are_almost_the_same(alpha_complex_from_points.filtration(f_simplex), 0.0));
        break;
      case 1:
        BOOST_CHECK(are_almost_the_same(alpha_complex_from_points.filtration(f_simplex), 1.0/2.0));
        break;
      case 2:
        BOOST_CHECK(are_almost_the_same(alpha_complex_from_points.filtration(f_simplex), 2.0/3.0));
        break;
      case 3:
        BOOST_CHECK(are_almost_the_same(alpha_complex_from_points.filtration(f_simplex), 3.0/4.0));
        break;
      default:
        BOOST_CHECK(false); // Shall not happen
        break;
    }
  }

}