PL_approximation.h
1 /* This file is part of the Gudhi Library - https://gudhi.inria.fr/ - which is released under MIT.
2  * See file LICENSE or go to https://gudhi.inria.fr/licensing/ for full license details.
3  * Author(s): Siargey Kachanovich
4  *
5  * Copyright (C) 2019 Inria
6  *
7  * Modification(s):
8  * - YYYY/MM Author: Description of the modification
9  */
10 
11 #ifndef FUNCTIONS_PL_APPROXIMATION_H_
12 #define FUNCTIONS_PL_APPROXIMATION_H_
13 
14 #include <cstdlib> // for std::size_t
15 
16 #include <Eigen/Dense>
17 
18 namespace Gudhi {
19 
20 namespace coxeter_triangulation {
21 
32 template <class Function_, class Triangulation_>
38  Eigen::VectorXd operator()(const Eigen::VectorXd& p) const {
39  std::size_t cod_d = this->cod_d();
40  std::size_t amb_d = this->amb_d();
41  auto s = tr_.locate_point(p);
42  Eigen::MatrixXd matrix(cod_d, s.dimension() + 1);
43  Eigen::MatrixXd vertex_matrix(amb_d + 1, s.dimension() + 1);
44  for (std::size_t i = 0; i < s.dimension() + 1; ++i) vertex_matrix(0, i) = 1;
45  std::size_t j = 0;
46  for (auto v : s.vertex_range()) {
47  Eigen::VectorXd pt_v = tr_.cartesian_coordinates(v);
48  Eigen::VectorXd fun_v = fun_(pt_v);
49  for (std::size_t i = 1; i < amb_d + 1; ++i) vertex_matrix(i, j) = pt_v(i - 1);
50  for (std::size_t i = 0; i < cod_d; ++i) matrix(i, j) = fun_v(i);
51  j++;
52  }
53  assert(j == s.dimension() + 1);
54  Eigen::VectorXd z(amb_d + 1);
55  z(0) = 1;
56  for (std::size_t i = 1; i < amb_d + 1; ++i) z(i) = p(i - 1);
57  Eigen::VectorXd lambda = vertex_matrix.colPivHouseholderQr().solve(z);
58  Eigen::VectorXd result = matrix * lambda;
59  return result;
60  }
61 
63  std::size_t amb_d() const { return fun_.amb_d(); }
64 
66  std::size_t cod_d() const { return fun_.cod_d(); }
67 
69  Eigen::VectorXd seed() const {
70  // TODO: not finished. Should use an oracle.
71  return Eigen::VectorXd(amb_d());
72  }
73 
81  PL_approximation(const Function_& function, const Triangulation_& triangulation)
82  : fun_(function), tr_(triangulation) {}
83 
84  private:
85  Function_ fun_;
86  Triangulation_ tr_;
87 };
88 
101 template <class Function_, class Triangulation_>
103  const Triangulation_& triangulation) {
104  return PL_approximation<Function_, Triangulation_>(function, triangulation);
105 }
106 
107 } // namespace coxeter_triangulation
108 
109 } // namespace Gudhi
110 
111 #endif
PL_approximation< Function_, Triangulation_ > make_pl_approximation(const Function_ &function, const Triangulation_ &triangulation)
Static constructor of the piecewise-linear approximation of a function induced by an ambient triangul...
Definition: PL_approximation.h:102
Gudhi namespace.
Definition: SimplicialComplexForAlpha.h:14
Constructs a piecewise-linear approximation of a function induced by an ambient triangulation.
Definition: PL_approximation.h:33
std::size_t amb_d() const
Returns the domain (ambient) dimension.
Definition: PL_approximation.h:63
Eigen::VectorXd seed() const
Returns a point on the zero-set.
Definition: PL_approximation.h:69
std::size_t cod_d() const
Returns the codomain dimension.
Definition: PL_approximation.h:66
Eigen::VectorXd operator()(const Eigen::VectorXd &p) const
Value of the function at a specified point.
Definition: PL_approximation.h:38
PL_approximation(const Function_ &function, const Triangulation_ &triangulation)
Constructor of the piecewise-linear approximation of a function induced by an ambient triangulation.
Definition: PL_approximation.h:81