summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorvrouvrea <vrouvrea@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-04-24 12:31:44 +0000
committervrouvrea <vrouvrea@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-04-24 12:31:44 +0000
commit3ce013afc21df5d05d663dd17a6640cec9af4aed (patch)
tree21c52342ecacc57b30be3aed653d3d4e9252d629
parentec96244695c80f6282c378a7d3d920677da4658d (diff)
Fix Python tangential complex version
git-svn-id: svn+ssh://scm.gforge.inria.fr/svnroot/gudhi/branches/python_2.1.0_fix_vincent@3392 636b058d-ea47-450e-bf9e-a15bfbe3eedb Former-commit-id: f0d9e008229a6cd7868f4afd11eaf3c566e81adf
-rw-r--r--src/cython/CMakeLists.txt2
-rw-r--r--src/cython/cython/tangential_complex.pyx15
-rw-r--r--src/cython/doc/tangential_complex_user.rst7
-rwxr-xr-xsrc/cython/example/tangential_complex_plain_homology_from_off_file_example.py5
-rw-r--r--src/cython/include/Tangential_complex_interface.h13
-rwxr-xr-xsrc/cython/test/test_tangential_complex.py2
6 files changed, 22 insertions, 22 deletions
diff --git a/src/cython/CMakeLists.txt b/src/cython/CMakeLists.txt
index b19cc550..f9721b46 100644
--- a/src/cython/CMakeLists.txt
+++ b/src/cython/CMakeLists.txt
@@ -183,7 +183,7 @@ if(CYTHON_FOUND)
WORKING_DIRECTORY ${CMAKE_CURRENT_BINARY_DIR}
COMMAND ${CMAKE_COMMAND} -E env "PYTHONPATH=${CMAKE_CURRENT_BINARY_DIR}"
${PYTHON_EXECUTABLE} "${CMAKE_CURRENT_SOURCE_DIR}/example/tangential_complex_plain_homology_from_off_file_example.py"
- --no-diagram -f ${CMAKE_SOURCE_DIR}/data/points/tore3D_300.off)
+ --no-diagram -i 2 -f ${CMAKE_SOURCE_DIR}/data/points/tore3D_300.off)
add_gudhi_py_test(test_tangential_complex)
diff --git a/src/cython/cython/tangential_complex.pyx b/src/cython/cython/tangential_complex.pyx
index d55bb050..44a3a96a 100644
--- a/src/cython/cython/tangential_complex.pyx
+++ b/src/cython/cython/tangential_complex.pyx
@@ -33,9 +33,9 @@ __license__ = "GPL v3"
cdef extern from "Tangential_complex_interface.h" namespace "Gudhi":
cdef cppclass Tangential_complex_interface "Gudhi::tangential_complex::Tangential_complex_interface":
- Tangential_complex_interface(vector[vector[double]] points)
+ Tangential_complex_interface(int intrisic_dim, vector[vector[double]] points)
# bool from_file is a workaround for cython to find the correct signature
- Tangential_complex_interface(string off_file, bool from_file)
+ Tangential_complex_interface(int intrisic_dim, string off_file, bool from_file)
vector[double] get_point(unsigned vertex)
unsigned number_of_vertices()
unsigned number_of_simplices()
@@ -54,9 +54,12 @@ cdef class TangentialComplex:
cdef Tangential_complex_interface * thisptr
# Fake constructor that does nothing but documenting the constructor
- def __init__(self, points=None, off_file=''):
+ def __init__(self, intrisic_dim, points=None, off_file=''):
"""TangentialComplex constructor.
+ :param intrisic_dim: Intrinsic dimension of the manifold.
+ :type intrisic_dim: integer
+
:param points: A list of points in d-Dimension.
:type points: list of list of double
@@ -67,17 +70,17 @@ cdef class TangentialComplex:
"""
# The real cython constructor
- def __cinit__(self, points=None, off_file=''):
+ def __cinit__(self, intrisic_dim, points=None, off_file=''):
if off_file is not '':
if os.path.isfile(off_file):
- self.thisptr = new Tangential_complex_interface(str.encode(off_file), True)
+ self.thisptr = new Tangential_complex_interface(intrisic_dim, str.encode(off_file), True)
else:
print("file " + off_file + " not found.")
else:
if points is None:
# Empty tangential construction
points=[]
- self.thisptr = new Tangential_complex_interface(points)
+ self.thisptr = new Tangential_complex_interface(intrisic_dim, points)
def __dealloc__(self):
diff --git a/src/cython/doc/tangential_complex_user.rst b/src/cython/doc/tangential_complex_user.rst
index efa6d7ce..fafb3193 100644
--- a/src/cython/doc/tangential_complex_user.rst
+++ b/src/cython/doc/tangential_complex_user.rst
@@ -122,8 +122,8 @@ This example builds the Tangential complex of point set read in an OFF file.
.. testcode::
import gudhi
- tc = gudhi.TangentialComplex(off_file=gudhi.__root_source_dir__ + \
- '/data/points/alphacomplexdoc.off')
+ tc = gudhi.TangentialComplex(intrisic_dim = 1,
+ off_file=gudhi.__root_source_dir__ + '/data/points/alphacomplexdoc.off')
result_str = 'Tangential contains ' + repr(tc.num_simplices()) + \
' simplices - ' + repr(tc.num_vertices()) + ' vertices.'
print(result_str)
@@ -169,7 +169,8 @@ simplices.
.. testcode::
import gudhi
- tc = gudhi.TangentialComplex(points=[[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]])
+ tc = gudhi.TangentialComplex(intrisic_dim = 1,
+ points=[[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]])
result_str = 'Tangential contains ' + repr(tc.num_vertices()) + ' vertices.'
print(result_str)
diff --git a/src/cython/example/tangential_complex_plain_homology_from_off_file_example.py b/src/cython/example/tangential_complex_plain_homology_from_off_file_example.py
index 6145e7f2..08ae5d07 100755
--- a/src/cython/example/tangential_complex_plain_homology_from_off_file_example.py
+++ b/src/cython/example/tangential_complex_plain_homology_from_off_file_example.py
@@ -33,10 +33,11 @@ parser = argparse.ArgumentParser(description='TangentialComplex creation from '
'points read in a OFF file.',
epilog='Example: '
'example/tangential_complex_plain_homology_from_off_file_example.py '
- '-f ../data/points/tore3D_300.off'
+ '-f ../data/points/tore3D_300.off -i 3'
'- Constructs a tangential complex with the '
'points from the given OFF file')
parser.add_argument("-f", "--file", type=str, required=True)
+parser.add_argument("-i", "--intrisic_dim", type=int, required=True)
parser.add_argument("-b", "--band_boot", type=float, default=0.)
parser.add_argument('--no-diagram', default=False, action='store_true' , help='Flag for not to display the diagrams')
@@ -48,7 +49,7 @@ with open(args.file, 'r') as f:
print("#####################################################################")
print("TangentialComplex creation from points read in a OFF file")
- tc = gudhi.TangentialComplex(off_file=args.file)
+ tc = gudhi.TangentialComplex(intrisic_dim = args.intrisic_dim, off_file=args.file)
st = tc.create_simplex_tree()
message = "Number of simplices=" + repr(st.num_simplices())
diff --git a/src/cython/include/Tangential_complex_interface.h b/src/cython/include/Tangential_complex_interface.h
index 0c3a510e..aef50424 100644
--- a/src/cython/include/Tangential_complex_interface.h
+++ b/src/cython/include/Tangential_complex_interface.h
@@ -45,24 +45,19 @@ class Tangential_complex_interface {
using TC = Tangential_complex<Dynamic_kernel, CGAL::Dynamic_dimension_tag, CGAL::Parallel_tag>;
public:
- Tangential_complex_interface(const std::vector<std::vector<double>>& points) {
+ Tangential_complex_interface(int intrisic_dim, const std::vector<std::vector<double>>& points) {
Dynamic_kernel k;
- unsigned intrisic_dim = 0;
- if (points.size() > 0)
- intrisic_dim = points[0].size() - 1;
tangential_complex_ = new TC(points, intrisic_dim, k);
tangential_complex_->compute_tangential_complex();
num_inconsistencies_ = tangential_complex_->number_of_inconsistent_simplices();
}
- Tangential_complex_interface(const std::string& off_file_name, bool from_file = true) {
- Gudhi::Points_off_reader<Point_d> off_reader(off_file_name);
+ Tangential_complex_interface(int intrisic_dim, const std::string& off_file_name, bool from_file = true) {
Dynamic_kernel k;
- unsigned intrisic_dim = 0;
+
+ Gudhi::Points_off_reader<Point_d> off_reader(off_file_name);
std::vector<Point_d> points = off_reader.get_point_cloud();
- if (points.size() > 0)
- intrisic_dim = points[0].size() - 1;
tangential_complex_ = new TC(points, intrisic_dim, k);
tangential_complex_->compute_tangential_complex();
diff --git a/src/cython/test/test_tangential_complex.py b/src/cython/test/test_tangential_complex.py
index 8aa4023c..ca9d936d 100755
--- a/src/cython/test/test_tangential_complex.py
+++ b/src/cython/test/test_tangential_complex.py
@@ -29,7 +29,7 @@ __license__ = "GPL v3"
def test_tangential():
point_list = [[0.0, 0.0], [1.0, 0.0], [0.0, 1.0], [1.0, 1.0]]
- tc = TangentialComplex(points=point_list)
+ tc = TangentialComplex(intrisic_dim = 1, points=point_list)
assert tc.__is_defined() == True
assert tc.num_vertices() == 4