summaryrefslogtreecommitdiff
path: root/src/Coxeter_triangulation/test/freud_triang_test.cpp
blob: 2cf8f00e40dd07b7191358e2af2160f7951d4369 (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
/*    This file is part of the Gudhi Library - https://gudhi.inria.fr/ - which is released under MIT.
 *    See file LICENSE or go to https://gudhi.inria.fr/licensing/ for full license details.
 *    Author(s):       Siargey Kachanovich
 *
 *    Copyright (C) 2019 Inria
 *
 *    Modification(s):
 *      - YYYY/MM Author: Description of the modification
 */

#define BOOST_TEST_DYN_LINK
#define BOOST_TEST_MODULE "freudenthal_triangulation"
#include <boost/test/unit_test.hpp>

#include <gudhi/Unitary_tests_utils.h>
#include <gudhi/Freudenthal_triangulation.h>
#include <gudhi/Coxeter_triangulation.h>

BOOST_AUTO_TEST_CASE(freudenthal_triangulation) {
  // Point location check
  typedef std::vector<double> Point;
  typedef Gudhi::coxeter_triangulation::Freudenthal_triangulation<> FK_triangulation;
  typedef typename FK_triangulation::Simplex_handle Simplex_handle;
  typedef typename FK_triangulation::Vertex_handle Vertex_handle;
  typedef typename Simplex_handle::OrderedSetPartition Ordered_set_partition;
  typedef typename Ordered_set_partition::value_type Part;

  FK_triangulation tr(3);

  // Point location check
  {
    Point point({3, -1, 0});
    Simplex_handle s = tr.locate_point(point);
    BOOST_CHECK(s.vertex() == Vertex_handle({3, -1, 0}));
    BOOST_CHECK(s.partition() == Ordered_set_partition({Part({0, 1, 2, 3})}));
  }

  {
    Point point({3.5, -1.5, 0.5});
    Simplex_handle s = tr.locate_point(point);
    BOOST_CHECK(s.vertex() == Vertex_handle({3, -2, 0}));
    BOOST_CHECK(s.partition() == Ordered_set_partition({Part({0, 1, 2}), Part({3})}));
  }

  {
    Point point({3.5, -1.8, 0.5});
    Simplex_handle s = tr.locate_point(point);
    BOOST_CHECK(s.vertex() == Vertex_handle({3, -2, 0}));
    BOOST_CHECK(s.partition() == Ordered_set_partition({Part({0, 2}), Part({1}), Part({3})}));
  }

  {
    Point point({3.5, -1.8, 0.3});
    Simplex_handle s = tr.locate_point(point);
    BOOST_CHECK(s.vertex() == Vertex_handle({3, -2, 0}));
    BOOST_CHECK(s.partition() == Ordered_set_partition({Part({0}), Part({2}), Part({1}), Part({3})}));
  }

  // Dimension check
  BOOST_CHECK(tr.dimension() == 3);
  // Matrix check
  Eigen::MatrixXd default_matrix = Eigen::MatrixXd::Identity(3, 3);
  BOOST_CHECK(tr.matrix() == default_matrix);
  // Vector check
  Eigen::MatrixXd default_offset = Eigen::VectorXd::Zero(3);
  BOOST_CHECK(tr.offset() == default_offset);

  // Barycenter check
  Point point({3.5, -1.8, 0.3});
  Simplex_handle s = tr.locate_point(point);
  Eigen::Vector3d barycenter_cart = Eigen::Vector3d::Zero();
  for (auto v : s.vertex_range())
    for (std::size_t i = 0; i < v.size(); i++) barycenter_cart(i) += v[i];
  barycenter_cart /= 4.;  // simplex is three-dimensional
  Eigen::Vector3d barycenter = tr.barycenter(s);
  for (std::size_t i = 0; (long int)i < barycenter.size(); i++)
    GUDHI_TEST_FLOAT_EQUALITY_CHECK(barycenter(i), barycenter_cart(i), 1e-7);

  // Barycenter check for twice the scale
  s = tr.locate_point(point, 2);
  barycenter_cart = Eigen::Vector3d::Zero();
  for (auto v : s.vertex_range())
    for (std::size_t i = 0; i < v.size(); i++) barycenter_cart(i) += v[i];
  barycenter_cart /= 3.;  // simplex is now a two-dimensional face
  barycenter_cart /= 2.;  // scale
  barycenter = tr.barycenter(s, 2);
  for (std::size_t i = 0; (long int)i < barycenter.size(); i++)
    GUDHI_TEST_FLOAT_EQUALITY_CHECK(barycenter(i), barycenter_cart(i), 1e-7);

  // Matrix and offset change check
  Eigen::MatrixXd new_matrix(3, 3);
  new_matrix << 1, 0, 0, -1, 1, 0, -1, 0, 1;
  Eigen::Vector3d new_offset(1.5, 1, 0.5);
  tr.change_matrix(new_matrix);
  tr.change_offset(new_offset);

  BOOST_CHECK(tr.matrix() == new_matrix);
  BOOST_CHECK(tr.offset() == new_offset);
}

#ifdef GUDHI_DEBUG
BOOST_AUTO_TEST_CASE(freudenthal_triangulation_exceptions_in_debug_mode) {
  // Point location check
  typedef Gudhi::coxeter_triangulation::Freudenthal_triangulation<> FK_triangulation;

  BOOST_CHECK_THROW (FK_triangulation tr(3, Eigen::MatrixXd::Identity(3, 3), Eigen::VectorXd::Zero(4)),
                     std::invalid_argument);

  FK_triangulation tr(3);
  // Point of dimension 4
  std::vector<double> point({3.5, -1.8, 0.3, 4.1});
  BOOST_CHECK_THROW (tr.locate_point(point), std::invalid_argument);
}
#endif