Loading...
Searching...
No Matches
IntersectionOracle.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 CONCEPT_COXETER_TRIANGULATION_INTERSECTION_ORACLE_H_
12#define CONCEPT_COXETER_TRIANGULATION_INTERSECTION_ORACLE_H_
13
14#include <cstdlib> // for std::size_t
15
16#include <Eigen/Dense>
17
18namespace Gudhi {
19
20namespace coxeter_triangulation {
21
28 std::size_t amb_d() const;
29
31 std::size_t cod_d() const;
32
53 template <class Simplex_handle, class Triangulation>
54 Query_result<Simplex_handle> intersects(const Simplex_handle& simplex, const Triangulation& triangulation) const;
55
76 template <class Simplex_handle, class Triangulation>
78 const Triangulation& triangulation) const;
79
90 template <class Triangulation>
91 bool lies_in_domain(const Eigen::VectorXd& p, const Triangulation& triangulation) const {
92 Eigen::VectorXd pl_p = make_pl_approximation(domain_fun_, triangulation)(p);
93 return pl_p(0) < 0;
94 }
95
97 const Function_& function() const;
98};
99
100} // namespace coxeter_triangulation
101
102} // namespace Gudhi
103
104#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
The concept IntersectionOracle describes the requirements for a type to implement an intersection ora...
Definition: IntersectionOracle.h:26
std::size_t cod_d() const
Returns the codomain dimension of the underlying manifold.
Query_result< Simplex_handle > intersects(const Simplex_handle &simplex, const Triangulation &triangulation) const
Intersection query with the relative interior of the manifold.
const Function_ & function() const
Returns the function that defines the interior of the manifold.
Query_result< Simplex_handle > intersects_boundary(const Simplex_handle &simplex, const Triangulation &triangulation) const
Intersection query with the boundary of the manifold.
bool lies_in_domain(const Eigen::VectorXd &p, const Triangulation &triangulation) const
Returns true if the input point lies inside the piecewise-linear domain induced by the given ambient ...
Definition: IntersectionOracle.h:91
std::size_t amb_d() const
Returns the domain (ambient) dimension of the underlying manifold.
The result of a query by an oracle such as Implicit_manifold_intersection_oracle.
Definition: Query_result.h:26