From 593131f9a55af9f82c7cfd8d272e1bb7eb0322e4 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Tue, 16 Aug 2016 20:10:21 +0200 Subject: fixed format selection for distance matrix --- ripser.cpp | 2 +- 1 file changed, 1 insertion(+), 1 deletion(-) diff --git a/ripser.cpp b/ripser.cpp index 1745927..168a724 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -861,7 +861,7 @@ int main(int argc, char** argv) { else if (parameter == "upper-distance") format = UPPER_DISTANCE_MATRIX; else if (parameter == "distance") - format = UPPER_DISTANCE_MATRIX; + format = DISTANCE_MATRIX; else if (parameter == "point-cloud") format = POINT_CLOUD; else if (parameter == "dipha") -- cgit v1.2.3 From 0a984f22f4c90c7a9d45b57187a6db128558d12c Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Tue, 16 Aug 2016 21:12:23 +0200 Subject: fixed csv reader --- ripser.cpp | 12 +++++++++--- 1 file changed, 9 insertions(+), 3 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index 168a724..6869aff 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -701,8 +701,11 @@ compressed_lower_distance_matrix read_point_cloud(std::istream& input_stream) { while (std::getline(input_stream, line)) { std::vector point; std::istringstream s(line); - while (s >> value) point.push_back(value); - if (!point.empty()) points.push_back(point); + while (s >> value) { + point.push_back(value); + s.ignore(); + } + if (!point.empty()) points.push_back(point); assert(point.size() == points.front().size()); } @@ -751,7 +754,10 @@ compressed_lower_distance_matrix read_distance_matrix(std::istream& input_stream value_t value; for (int i = 0; std::getline(input_stream, line); ++i) { std::istringstream s(line); - for (int j = 0; j < i && s >> value; ++j) distances.push_back(value); + for (int j = 0; j < i && s >> value; ++j) { + distances.push_back(value); + s.ignore(); + } } return compressed_lower_distance_matrix(std::move(distances)); -- cgit v1.2.3 From c130a8b6dd389456eb7d09f53f6fe676358db3d6 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 18 Aug 2016 18:21:08 +0200 Subject: added reduction matrix option to Makefile updated readme --- Makefile | 7 +++++-- README.md | 4 ++-- 2 files changed, 7 insertions(+), 4 deletions(-) diff --git a/Makefile b/Makefile index 9359d43..331fdd6 100644 --- a/Makefile +++ b/Makefile @@ -1,7 +1,7 @@ build: ripser -all: ripser ripser-coeff +all: ripser ripser-coeff ripser-reduction ripser: ripser.cpp @@ -10,6 +10,9 @@ ripser: ripser.cpp ripser-coeff: ripser.cpp c++ -std=c++11 ripser.cpp -o ripser-coeff -Ofast -D NDEBUG -D USE_COEFFICIENTS +ripser-reduction: ripser.cpp + c++ -std=c++11 ripser.cpp -o ripser-reduction -Ofast -D NDEBUG -D ASSEMBLE_REDUCTION_MATRIX + clean: - rm ripser ripser-coeff \ No newline at end of file + rm -f ripser ripser-coeff ripser-reduction diff --git a/README.md b/README.md index 13ae86b..1f41a15 100644 --- a/README.md +++ b/README.md @@ -53,7 +53,7 @@ make Ripser supports several compile-time options. They are switched on by defining the C preprocessor macros listed below, either using `#define` in the code or by passing an argument to the compiler. The following options are supported: - - `ASSEMBLE_REDUCTION_MATRIX`: store the reduction matrix; may speed up computation but will also increase memory usage + - `ASSEMBLE_REDUCTION_MATRIX`: store the reduction matrix; may affect computation time but also memory usage; recommended for large and difficult problem instances - `USE_COEFFICIENTS`: enable support for coefficients in a prime field - `INDICATE_PROGRESS`: indicate the current progress in the console - `PRINT_PERSISTENCE_PAIRS`: output the computed persistence pairs (enabled by default in the code; comment out to disable) @@ -65,7 +65,7 @@ For example, to build Ripser with support for coefficients: $ c++ -std=c++11 ripser.cpp -o ripser -Ofast -D NDEBUG -D USE_COEFFICIENTS ``` -A Makefile is provided with some variants of the above options. Use `make all` to build them. The default `make` builds a binary without any of the above option. +A Makefile is provided with some variants of the above options. Use `make all` to build them. The default `make` builds a binary with the default options. The input is given either in a file whose name is passed as an argument, or through stdin. The following options are supported at the command line: -- cgit v1.2.3 From cb53a8a3108735a4ba91afb8a8084bccf1d65f06 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Sun, 21 Aug 2016 16:24:07 +0200 Subject: cleanup to remove compiler warning --- ripser.cpp | 19 +++++++++---------- 1 file changed, 9 insertions(+), 10 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index 6869aff..200ada1 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -251,14 +251,13 @@ public: class simplex_coboundary_enumerator { private: - index_t idx, modified_idx, dim, v, k; + index_t idx, modified_idx, v, k; const binomial_coeff_table& binomial_coeff; public: simplex_coboundary_enumerator(index_t _idx, index_t _dim, index_t _n, const binomial_coeff_table& _binomial_coeff) - : idx(_idx), modified_idx(_idx), dim(_dim), k(dim + 1), v(_n - 1), - binomial_coeff(_binomial_coeff) {} + : idx(_idx), modified_idx(_idx), v(_n - 1), k(_dim + 1), binomial_coeff(_binomial_coeff) {} bool has_next() { while ((v != -1) && (binomial_coeff(v, k) <= idx)) { @@ -702,10 +701,10 @@ compressed_lower_distance_matrix read_point_cloud(std::istream& input_stream) { std::vector point; std::istringstream s(line); while (s >> value) { - point.push_back(value); - s.ignore(); - } - if (!point.empty()) points.push_back(point); + point.push_back(value); + s.ignore(); + } + if (!point.empty()) points.push_back(point); assert(point.size() == points.front().size()); } @@ -755,9 +754,9 @@ compressed_lower_distance_matrix read_distance_matrix(std::istream& input_stream for (int i = 0; std::getline(input_stream, line); ++i) { std::istringstream s(line); for (int j = 0; j < i && s >> value; ++j) { - distances.push_back(value); - s.ignore(); - } + distances.push_back(value); + s.ignore(); + } } return compressed_lower_distance_matrix(std::move(distances)); -- cgit v1.2.3 From a76f7e89dd7f9833265147964d2f10c4cfd7e7bf Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Sun, 21 Aug 2016 22:21:06 +0200 Subject: cleanup --- .clang-format | 2 +- ripser.cpp | 131 +++++++++++++++++++++------------------------------------- 2 files changed, 47 insertions(+), 86 deletions(-) diff --git a/.clang-format b/.clang-format index 1975b19..598d7c2 100644 --- a/.clang-format +++ b/.clang-format @@ -1,7 +1,7 @@ BasedOnStyle: LLVM IndentWidth: 4 TabWidth: 4 -ColumnLimit: 100 +ColumnLimit: 120 AccessModifierOffset: -4 AllowShortIfStatementsOnASingleLine: true AllowShortLoopsOnASingleLine: true diff --git a/ripser.cpp b/ripser.cpp index 200ada1..d3c5127 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -96,8 +96,7 @@ std::vector multiplicative_inverse_vector(const coefficient_t m) template OutputIterator get_simplex_vertices(index_t idx, const index_t dim, index_t n, - const binomial_coeff_table& binomial_coeff, - OutputIterator out) { + const binomial_coeff_table& binomial_coeff, OutputIterator out) { --n; for (index_t k = dim + 1; k > 0; --k) { @@ -124,8 +123,7 @@ OutputIterator get_simplex_vertices(index_t idx, const index_t dim, index_t n, return out; } -std::vector vertices_of_simplex(const index_t simplex_index, const index_t dim, - const index_t n, +std::vector vertices_of_simplex(const index_t simplex_index, const index_t dim, const index_t n, const binomial_coeff_table& binomial_coeff) { std::vector vertices; get_simplex_vertices(simplex_index, dim, n, binomial_coeff, std::back_inserter(vertices)); @@ -136,15 +134,12 @@ std::vector vertices_of_simplex(const index_t simplex_index, const inde struct entry_t { index_t index; coefficient_t coefficient; - entry_t(index_t _index, coefficient_t _coefficient) - : index(_index), coefficient(_coefficient) {} + entry_t(index_t _index, coefficient_t _coefficient) : index(_index), coefficient(_coefficient) {} entry_t(index_t _index) : index(_index), coefficient(1) {} entry_t() : index(0), coefficient(1) {} }; -entry_t make_entry(index_t _index, coefficient_t _coefficient) { - return entry_t(_index, _coefficient); -} +entry_t make_entry(index_t _index, coefficient_t _coefficient) { return entry_t(_index, _coefficient); } index_t get_index(entry_t e) { return e.index; } index_t get_coefficient(entry_t e) { return e.coefficient; } void set_coefficient(entry_t& e, const coefficient_t c) { e.coefficient = c; } @@ -188,20 +183,14 @@ public: const entry_t& get_entry(const diameter_entry_t& p) { return p.second; } entry_t& get_entry(diameter_entry_t& p) { return p.second; } const index_t get_index(const diameter_entry_t& p) { return get_index(get_entry(p)); } -const coefficient_t get_coefficient(const diameter_entry_t& p) { - return get_coefficient(get_entry(p)); -} +const coefficient_t get_coefficient(const diameter_entry_t& p) { return get_coefficient(get_entry(p)); } const value_t& get_diameter(const diameter_entry_t& p) { return p.first; } -void set_coefficient(diameter_entry_t& p, const coefficient_t c) { - set_coefficient(get_entry(p), c); -} -diameter_entry_t make_diameter_entry(value_t _diameter, index_t _index, - coefficient_t _coefficient) { +void set_coefficient(diameter_entry_t& p, const coefficient_t c) { set_coefficient(get_entry(p), c); } +diameter_entry_t make_diameter_entry(value_t _diameter, index_t _index, coefficient_t _coefficient) { return std::make_pair(_diameter, make_entry(_index, _coefficient)); } diameter_entry_t make_diameter_entry(diameter_index_t _diameter_index, coefficient_t _coefficient) { - return std::make_pair(get_diameter(_diameter_index), - make_entry(get_index(_diameter_index), _coefficient)); + return std::make_pair(get_diameter(_diameter_index), make_entry(get_index(_diameter_index), _coefficient)); } template struct greater_diameter_or_smaller_index { @@ -230,9 +219,7 @@ public: get_simplex_vertices(index, dim, dist.size(), binomial_coeff, vertices.begin()); for (index_t i = 0; i <= dim; ++i) - for (index_t j = 0; j < i; ++j) { - diam = std::max(diam, dist(vertices[i], vertices[j])); - } + for (index_t j = 0; j < i; ++j) { diam = std::max(diam, dist(vertices[i], vertices[j])); } return diam; } @@ -240,8 +227,8 @@ public: assert(a < binomial_coeff(dist.size(), dim + 1)); assert(b < binomial_coeff(dist.size(), dim + 1)); - return greater_diameter_or_smaller_index()( - diameter_index_t(diameter(a), a), diameter_index_t(diameter(b), b)); + return greater_diameter_or_smaller_index()(diameter_index_t(diameter(a), a), + diameter_index_t(diameter(b), b)); } template bool operator()(const Entry& a, const Entry& b) const { @@ -255,8 +242,7 @@ private: const binomial_coeff_table& binomial_coeff; public: - simplex_coboundary_enumerator(index_t _idx, index_t _dim, index_t _n, - const binomial_coeff_table& _binomial_coeff) + simplex_coboundary_enumerator(index_t _idx, index_t _dim, index_t _n, const binomial_coeff_table& _binomial_coeff) : idx(_idx), modified_idx(_idx), v(_n - 1), k(_dim + 1), binomial_coeff(_binomial_coeff) {} bool has_next() { @@ -271,8 +257,7 @@ public: } std::pair next() { - auto result = - std::make_pair(make_entry(modified_idx + binomial_coeff(v, k + 1), k & 1 ? -1 : 1), v); + auto result = std::make_pair(make_entry(modified_idx + binomial_coeff(v, k + 1), k & 1 ? -1 : 1), v); --v; return result; } @@ -323,16 +308,14 @@ template <> void compressed_distance_matrix::init_rows() { } } -template <> -value_t compressed_distance_matrix::operator()(const index_t i, - const index_t j) const { - return i == j ? 0 : rows[std::min(i, j)][std::max(i, j)]; +template <> value_t compressed_distance_matrix::operator()(index_t i, index_t j) const { + std::tie(i, j) = std::minmax(i, j); + return i == j ? 0 : rows[i][j]; } -template <> -value_t compressed_distance_matrix::operator()(const index_t i, - const index_t j) const { - return i == j ? 0 : rows[std::max(i, j)][std::min(i, j)]; +template <> value_t compressed_distance_matrix::operator()(index_t i, index_t j) const { + std::tie(i, j) = std::minmax(i, j); + return i == j ? 0 : rows[j][i]; } typedef compressed_distance_matrix compressed_lower_distance_matrix; @@ -345,9 +328,9 @@ public: euclidean_distance_matrix(std::vector>&& _points) : points(_points) {} value_t operator()(const index_t i, const index_t j) const { - return std::sqrt(std::inner_product( - points[i].begin(), points[i].end(), points[j].begin(), value_t(), std::plus(), - [](value_t u, value_t v) { return (u - v) * (u - v); })); + return std::sqrt(std::inner_product(points[i].begin(), points[i].end(), points[j].begin(), value_t(), + std::plus(), + [](value_t u, value_t v) { return (u - v) * (u - v); })); } size_t size() const { return points.size(); } @@ -436,17 +419,15 @@ public: } }; -template -void push_entry(Heap& column, index_t i, coefficient_t c, value_t diameter) { +template void push_entry(Heap& column, index_t i, coefficient_t c, value_t diameter) { entry_t e = make_entry(i, c); column.push(std::make_pair(diameter, e)); } template void assemble_columns_to_reduce(std::vector& columns_to_reduce, - hash_map& pivot_column_index, - const Comparator& comp, index_t dim, index_t n, value_t threshold, - const binomial_coeff_table& binomial_coeff) { + hash_map& pivot_column_index, const Comparator& comp, index_t dim, + index_t n, value_t threshold, const binomial_coeff_table& binomial_coeff) { index_t num_simplices = binomial_coeff(n, dim + 2); columns_to_reduce.clear(); @@ -476,9 +457,8 @@ void assemble_columns_to_reduce(std::vector& columns_to_reduce } template -void compute_pairs(std::vector& columns_to_reduce, - hash_map& pivot_column_index, const DistanceMatrix& dist, - const ComparatorCofaces& comp, const Comparator& comp_prev, index_t dim, +void compute_pairs(std::vector& columns_to_reduce, hash_map& pivot_column_index, + const DistanceMatrix& dist, const ComparatorCofaces& comp, const Comparator& comp_prev, index_t dim, index_t n, value_t threshold, coefficient_t modulus, const std::vector& multiplicative_inverse, const binomial_coeff_table& binomial_coeff) { @@ -502,8 +482,7 @@ void compute_pairs(std::vector& columns_to_reduce, auto column_to_reduce = columns_to_reduce[i]; #ifdef ASSEMBLE_REDUCTION_MATRIX - std::priority_queue, - smaller_index> + std::priority_queue, smaller_index> reduction_column; #endif @@ -516,8 +495,8 @@ void compute_pairs(std::vector& columns_to_reduce, #ifdef INDICATE_PROGRESS if ((i + 1) % 1000 == 0) std::cout << "\033[K" - << "reducing column " << i + 1 << "/" << columns_to_reduce.size() - << " (diameter " << diameter << ")" << std::flush << "\r"; + << "reducing column " << i + 1 << "/" << columns_to_reduce.size() << " (diameter " << diameter + << ")" << std::flush << "\r"; #endif index_t j = i; @@ -565,13 +544,11 @@ void compute_pairs(std::vector& columns_to_reduce, assert(simplex_diameter == comp_prev.diameter(get_index(simplex))); #ifdef ASSEMBLE_REDUCTION_MATRIX - reduction_column.push( - make_diameter_entry(simplex_diameter, get_index(simplex), simplex_coefficient)); + reduction_column.push(make_diameter_entry(simplex_diameter, get_index(simplex), simplex_coefficient)); #endif vertices.clear(); - get_simplex_vertices(get_index(simplex), dim, n, binomial_coeff, - std::back_inserter(vertices)); + get_simplex_vertices(get_index(simplex), dim, n, binomial_coeff, std::back_inserter(vertices)); coface_entries.clear(); simplex_coboundary_enumerator cofaces(get_index(simplex), dim, n, binomial_coeff); @@ -581,14 +558,11 @@ void compute_pairs(std::vector& columns_to_reduce, index_t covertex = coface_descriptor.second; index_t coface_index = get_index(coface); value_t coface_diameter = simplex_diameter; - for (index_t v : vertices) { - coface_diameter = std::max(coface_diameter, dist(v, covertex)); - } + for (index_t v : vertices) { coface_diameter = std::max(coface_diameter, dist(v, covertex)); } assert(comp.diameter(coface_index) == coface_diameter); if (coface_diameter <= threshold) { - coefficient_t coface_coefficient = - get_coefficient(coface) * get_coefficient(simplex) * factor; + coefficient_t coface_coefficient = get_coefficient(coface) * get_coefficient(simplex) * factor; coface_coefficient %= modulus; if (coface_coefficient < 0) coface_coefficient += modulus; assert(coface_coefficient >= 0); @@ -660,8 +634,7 @@ void compute_pairs(std::vector& columns_to_reduce, #else const coefficient_t coefficient = 1; #endif - reduction_matrix.push_back( - make_diameter_entry(get_diameter(e), index, coefficient)); + reduction_matrix.push_back(make_diameter_entry(get_diameter(e), index, coefficient)); } #else #ifdef USE_COEFFICIENTS @@ -678,13 +651,7 @@ void compute_pairs(std::vector& columns_to_reduce, #endif } -enum file_format { - LOWER_DISTANCE_MATRIX, - UPPER_DISTANCE_MATRIX, - DISTANCE_MATRIX, - POINT_CLOUD, - DIPHA -}; +enum file_format { LOWER_DISTANCE_MATRIX, UPPER_DISTANCE_MATRIX, DISTANCE_MATRIX, POINT_CLOUD, DIPHA }; template T read(std::istream& s) { T result; @@ -712,8 +679,7 @@ compressed_lower_distance_matrix read_point_cloud(std::istream& input_stream) { index_t n = eucl_dist.size(); - std::cout << "point cloud with " << n << " points in dimension " - << eucl_dist.points.front().size() << std::endl; + std::cout << "point cloud with " << n << " points in dimension " << eucl_dist.points.front().size() << std::endl; std::vector distances; @@ -810,16 +776,12 @@ void print_usage_and_exit(int exit_code) { << "Options:" << std::endl << std::endl << " --help print this screen" << std::endl - << " --format use the specified file format for the input. Options are:" - << std::endl - << " lower-distance (lower triangular distance matrix; default)" - << std::endl - << " upper-distance (upper triangular distance matrix)" - << std::endl + << " --format use the specified file format for the input. Options are:" << std::endl + << " lower-distance (lower triangular distance matrix; default)" << std::endl + << " upper-distance (upper triangular distance matrix)" << std::endl << " distance (full distance matrix)" << std::endl << " point-cloud (point cloud in Euclidean space)" << std::endl - << " dipha (distance matrix in DIPHA file format)" - << std::endl + << " dipha (distance matrix in DIPHA file format)" << std::endl << " --dim compute persistent homology up to dimension " << std::endl << " --threshold compute Rips complexes up to diameter " << std::endl #ifdef USE_COEFFICIENTS @@ -898,8 +860,8 @@ int main(int argc, char** argv) { std::cout << "distance matrix with " << n << " points" << std::endl; - auto result = std::minmax_element(dist.distances.begin(), dist.distances.end()); - std::cout << "value range: [" << *result.first << "," << *result.second << "]" << std::endl; + auto value_range = std::minmax_element(dist.distances.begin(), dist.distances.end()); + std::cout << "value range: [" << *value_range.first << "," << *value_range.second << "]" << std::endl; dim_max = std::min(dim_max, n - 2); @@ -917,11 +879,10 @@ int main(int argc, char** argv) { hash_map pivot_column_index; pivot_column_index.reserve(columns_to_reduce.size()); - compute_pairs(columns_to_reduce, pivot_column_index, dist, comp, comp_prev, dim, n, - threshold, modulus, multiplicative_inverse, binomial_coeff); + compute_pairs(columns_to_reduce, pivot_column_index, dist, comp, comp_prev, dim, n, threshold, modulus, + multiplicative_inverse, binomial_coeff); if (dim < dim_max) - assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, - threshold, binomial_coeff); + assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, threshold, binomial_coeff); } } \ No newline at end of file -- cgit v1.2.3 From 97a760c54be95cb949d1e9c25f9e4513982f8109 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Wed, 24 Aug 2016 10:52:21 +0200 Subject: fixed compressed_distance_matrix access --- ripser.cpp | 4 ++-- 1 file changed, 2 insertions(+), 2 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index d3c5127..ec58cb2 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -309,12 +309,12 @@ template <> void compressed_distance_matrix::init_rows() { } template <> value_t compressed_distance_matrix::operator()(index_t i, index_t j) const { - std::tie(i, j) = std::minmax(i, j); + if (i > j) std::swap(i, j); return i == j ? 0 : rows[i][j]; } template <> value_t compressed_distance_matrix::operator()(index_t i, index_t j) const { - std::tie(i, j) = std::minmax(i, j); + if (i > j) std::swap(i, j); return i == j ? 0 : rows[j][i]; } -- cgit v1.2.3 From 4c0ccbe423f29feadc9f40645fe326c4f56a8920 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Wed, 24 Aug 2016 23:06:41 +0200 Subject: union find for H^0 --- ripser.cpp | 79 ++++++++++++++++++++++++++++++++++++++++++++++++++++++++++---- 1 file changed, 75 insertions(+), 4 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index ec58cb2..a80a985 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -33,6 +33,19 @@ #include #include +#include +#include + +#include "prettyprint.hpp" + + +#ifdef __EMSCRIPTEN__ +#include +#include + +using namespace emscripten; +#endif + #ifdef USE_GOOGLE_HASHMAP #include template class hash_map : public google::sparse_hash_map { @@ -559,7 +572,7 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map multiplicative_inverse(multiplicative_inverse_vector(modulus)); std::vector columns_to_reduce; - for (index_t index = n; index-- > 0;) columns_to_reduce.push_back(diameter_index_t(0, index)); + + boost::vector_property_map rank(n); + boost::vector_property_map parent(n); + + boost::disjoint_sets, boost::vector_property_map> dset(rank, parent); + + for (index_t index = n; index-- > 0;) { + dset.make_set(index); + } + + std::vector edges(binomial_coeff(n, 2)); + + rips_filtration_comparator comp(dist, 1, binomial_coeff); + + assert(dist.distances.size() == binomial_coeff(n, 2)); + + for (index_t index = binomial_coeff(n, 2); index-- > 0;) { + edges[index] = diameter_index_t(comp.diameter(index), index); + } + + std::sort(edges.rbegin(), edges.rend(), comp); + + std::vector vertices_of_edge(2); + +#ifdef PRINT_PERSISTENCE_PAIRS + std::cout << "persistence intervals in dim 0:" << std::endl; +#endif + + for (auto e: edges) { + + vertices_of_edge.clear(); + get_simplex_vertices(get_index(e), 1, n, binomial_coeff, std::back_inserter(vertices_of_edge)); + index_t u = dset.find_set(vertices_of_edge[0]), v = dset.find_set(vertices_of_edge[1]); + + if ( u != v ) { +// std::cout << "tree edge " << vertices_of_edge[0] << "-" << vertices_of_edge[1] << " " << e << std::endl; +#ifdef PRINT_PERSISTENCE_PAIRS +#ifdef INDICATE_PROGRESS + std::cout << "\033[K"; +#endif + std::cout << " [0," << get_diameter(e) << ")" << std::endl; +#endif + dset.link(u, v); + } else { +// std::cout << "cycle edge " << vertices_of_edge[0] << "-" << vertices_of_edge[1] << " " << e << std::endl; + columns_to_reduce.push_back(e); + } + } + + std::sort(columns_to_reduce.begin(), columns_to_reduce.end(), + greater_diameter_or_smaller_index()); + +// std::cout << columns_to_reduce << std::endl; + - for (index_t dim = 0; dim <= dim_max; ++dim) { +// columns_to_reduce.clear(); +// for (index_t index = n; index-- > 0;) columns_to_reduce.push_back(diameter_index_t(0, index)); + + for (index_t dim = 1; dim <= dim_max; ++dim) { rips_filtration_comparator comp(dist, dim + 1, binomial_coeff); rips_filtration_comparator comp_prev(dist, dim, binomial_coeff); @@ -882,7 +951,9 @@ int main(int argc, char** argv) { compute_pairs(columns_to_reduce, pivot_column_index, dist, comp, comp_prev, dim, n, threshold, modulus, multiplicative_inverse, binomial_coeff); - if (dim < dim_max) + if (dim < dim_max) { assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, threshold, binomial_coeff); +// std::cout << columns_to_reduce << std::endl; + } } } \ No newline at end of file -- cgit v1.2.3 From bedbfdfb53779c85d17324efbb2867d308aa174b Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Wed, 24 Aug 2016 23:51:52 +0200 Subject: union find implementation --- ripser.cpp | 56 +++++++++++++++++++++++++++++++++++++++++++++----------- 1 file changed, 45 insertions(+), 11 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index a80a985..fcf172c 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -33,8 +33,7 @@ #include #include -#include -#include + #include "prettyprint.hpp" @@ -62,6 +61,7 @@ typedef float value_t; typedef long index_t; typedef long coefficient_t; + class binomial_coeff_table { std::vector> B; index_t n_max, k_max; @@ -349,6 +349,47 @@ public: size_t size() const { return points.size(); } }; +class union_find +{ + std::vector parent; + std::vector rank; +public: + + union_find(index_t n): rank(n, 0), parent(n) + { + for (index_t i = 0; i < n; ++i) + parent[i] = i; + } + + index_t find(index_t x) { + index_t y = x, z = parent[y]; + while (z != y) { + y = z; + z = parent[y]; + } + y = parent[x]; + while (z != y) { + parent[x] = z; + x = y; + y = parent[x]; + } + return z; + } + void link(index_t x, index_t y) + { + x = find(x); + y = find(y); + if (x == y) return; + if (rank[x] > rank[y]) + parent[y] = x; + else { + parent[x] = y; + if (rank[x] == rank[y]) + ++rank[y]; + } + } +}; + template diameter_entry_t pop_pivot(Heap& column, coefficient_t modulus) { if (column.empty()) return diameter_entry_t(-1); @@ -883,14 +924,7 @@ int main(int argc, char** argv) { std::vector columns_to_reduce; - boost::vector_property_map rank(n); - boost::vector_property_map parent(n); - - boost::disjoint_sets, boost::vector_property_map> dset(rank, parent); - - for (index_t index = n; index-- > 0;) { - dset.make_set(index); - } + union_find dset(n); std::vector edges(binomial_coeff(n, 2)); @@ -914,7 +948,7 @@ int main(int argc, char** argv) { vertices_of_edge.clear(); get_simplex_vertices(get_index(e), 1, n, binomial_coeff, std::back_inserter(vertices_of_edge)); - index_t u = dset.find_set(vertices_of_edge[0]), v = dset.find_set(vertices_of_edge[1]); + index_t u = dset.find(vertices_of_edge[0]), v = dset.find(vertices_of_edge[1]); if ( u != v ) { // std::cout << "tree edge " << vertices_of_edge[0] << "-" << vertices_of_edge[1] << " " << e << std::endl; -- cgit v1.2.3 From 48b815399d656ea9988ede177c955000c8a67acf Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 25 Aug 2016 00:23:31 +0200 Subject: cleanup --- ripser.cpp | 97 +++++++++++++++++++++++--------------------------------------- 1 file changed, 35 insertions(+), 62 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index fcf172c..7a06590 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -33,18 +33,6 @@ #include #include - - -#include "prettyprint.hpp" - - -#ifdef __EMSCRIPTEN__ -#include -#include - -using namespace emscripten; -#endif - #ifdef USE_GOOGLE_HASHMAP #include template class hash_map : public google::sparse_hash_map { @@ -61,7 +49,6 @@ typedef float value_t; typedef long index_t; typedef long coefficient_t; - class binomial_coeff_table { std::vector> B; index_t n_max, k_max; @@ -613,7 +600,7 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map columns_to_reduce; - union_find dset(n); - - std::vector edges(binomial_coeff(n, 2)); - - rips_filtration_comparator comp(dist, 1, binomial_coeff); - - assert(dist.distances.size() == binomial_coeff(n, 2)); - - for (index_t index = binomial_coeff(n, 2); index-- > 0;) { - edges[index] = diameter_index_t(comp.diameter(index), index); - } - - std::sort(edges.rbegin(), edges.rend(), comp); - - std::vector vertices_of_edge(2); - + { + union_find dset(n); + std::vector edges; + rips_filtration_comparator comp(dist, 1, binomial_coeff); + for (index_t index = binomial_coeff(n, 2); index-- > 0;) { + value_t diameter = comp.diameter(index); + if (diameter <= threshold) edges.push_back(diameter_index_t(diameter, index)); + } + std::sort(edges.rbegin(), edges.rend(), greater_diameter_or_smaller_index()); + #ifdef PRINT_PERSISTENCE_PAIRS - std::cout << "persistence intervals in dim 0:" << std::endl; + std::cout << "persistence intervals in dim 0:" << std::endl; #endif - - for (auto e: edges) { - vertices_of_edge.clear(); - get_simplex_vertices(get_index(e), 1, n, binomial_coeff, std::back_inserter(vertices_of_edge)); - index_t u = dset.find(vertices_of_edge[0]), v = dset.find(vertices_of_edge[1]); - - if ( u != v ) { -// std::cout << "tree edge " << vertices_of_edge[0] << "-" << vertices_of_edge[1] << " " << e << std::endl; + std::vector vertices_of_edge(2); + for (auto e: edges) { + vertices_of_edge.clear(); + get_simplex_vertices(get_index(e), 1, n, binomial_coeff, std::back_inserter(vertices_of_edge)); + index_t u = dset.find(vertices_of_edge[0]), v = dset.find(vertices_of_edge[1]); + + if ( u != v ) { #ifdef PRINT_PERSISTENCE_PAIRS -#ifdef INDICATE_PROGRESS - std::cout << "\033[K"; -#endif - std::cout << " [0," << get_diameter(e) << ")" << std::endl; + std::cout << " [0," << get_diameter(e) << ")" << std::endl; #endif - dset.link(u, v); - } else { -// std::cout << "cycle edge " << vertices_of_edge[0] << "-" << vertices_of_edge[1] << " " << e << std::endl; - columns_to_reduce.push_back(e); + dset.link(u, v); + } else + columns_to_reduce.push_back(e); } + std::reverse(columns_to_reduce.begin(), columns_to_reduce.end()); + +#ifdef PRINT_PERSISTENCE_PAIRS + for (index_t i = 0; i < n; ++i) + if (dset.find(i) == i) + std::cout << " [0, )" << std::endl << std::flush; +#endif } - std::sort(columns_to_reduce.begin(), columns_to_reduce.end(), - greater_diameter_or_smaller_index()); - -// std::cout << columns_to_reduce << std::endl; - - - -// columns_to_reduce.clear(); -// for (index_t index = n; index-- > 0;) columns_to_reduce.push_back(diameter_index_t(0, index)); - for (index_t dim = 1; dim <= dim_max; ++dim) { rips_filtration_comparator comp(dist, dim + 1, binomial_coeff); rips_filtration_comparator comp_prev(dist, dim, binomial_coeff); - + hash_map pivot_column_index; pivot_column_index.reserve(columns_to_reduce.size()); - + compute_pairs(columns_to_reduce, pivot_column_index, dist, comp, comp_prev, dim, n, threshold, modulus, - multiplicative_inverse, binomial_coeff); - + multiplicative_inverse, binomial_coeff); + if (dim < dim_max) { assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, threshold, binomial_coeff); -// std::cout << columns_to_reduce << std::endl; + // std::cout << columns_to_reduce << std::endl; } } } \ No newline at end of file -- cgit v1.2.3 From c10327ed6ccf9e3c8961833cbf95431433d6b282 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 25 Aug 2016 00:36:02 +0200 Subject: cleanup --- ripser.cpp | 55 +++++++++++++++++++++++-------------------------------- 1 file changed, 23 insertions(+), 32 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index 7a06590..b4b528c 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -336,18 +336,15 @@ public: size_t size() const { return points.size(); } }; -class union_find -{ +class union_find { std::vector parent; std::vector rank; + public: - - union_find(index_t n): rank(n, 0), parent(n) - { - for (index_t i = 0; i < n; ++i) - parent[i] = i; + union_find(index_t n) : parent(n), rank(n, 0) { + for (index_t i = 0; i < n; ++i) parent[i] = i; } - + index_t find(index_t x) { index_t y = x, z = parent[y]; while (z != y) { @@ -362,8 +359,7 @@ public: } return z; } - void link(index_t x, index_t y) - { + void link(index_t x, index_t y) { x = find(x); y = find(y); if (x == y) return; @@ -371,8 +367,7 @@ public: parent[y] = x; else { parent[x] = y; - if (rank[x] == rank[y]) - ++rank[y]; + if (rank[x] == rank[y]) ++rank[y]; } } }; @@ -420,10 +415,10 @@ template diameter_entry_t get_pivot(Heap& column, coefficient_t } template class compressed_sparse_matrix { -public: std::vector bounds; std::vector entries; +public: size_t size() const { return bounds.size(); } typename std::vector::const_iterator cbegin(size_t index) const { @@ -528,8 +523,7 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map, - greater_diameter_or_smaller_index> - working_coboundary; + greater_diameter_or_smaller_index> working_coboundary; value_t diameter = get_diameter(column_to_reduce); @@ -850,9 +844,7 @@ int main(int argc, char** argv) { for (index_t i = 1; i < argc; ++i) { const std::string arg(argv[i]); - if (arg == "--help") { - print_usage_and_exit(0); - } else if (arg == "--dim") { + if (arg == "--help") { print_usage_and_exit(0); } else if (arg == "--dim") { std::string parameter = std::string(argv[++i]); size_t next_pos; dim_max = std::stol(parameter, &next_pos); @@ -910,7 +902,7 @@ int main(int argc, char** argv) { std::vector multiplicative_inverse(multiplicative_inverse_vector(modulus)); std::vector columns_to_reduce; - + { union_find dset(n); std::vector edges; @@ -920,18 +912,18 @@ int main(int argc, char** argv) { if (diameter <= threshold) edges.push_back(diameter_index_t(diameter, index)); } std::sort(edges.rbegin(), edges.rend(), greater_diameter_or_smaller_index()); - + #ifdef PRINT_PERSISTENCE_PAIRS std::cout << "persistence intervals in dim 0:" << std::endl; #endif - + std::vector vertices_of_edge(2); - for (auto e: edges) { + for (auto e : edges) { vertices_of_edge.clear(); get_simplex_vertices(get_index(e), 1, n, binomial_coeff, std::back_inserter(vertices_of_edge)); index_t u = dset.find(vertices_of_edge[0]), v = dset.find(vertices_of_edge[1]); - - if ( u != v ) { + + if (u != v) { #ifdef PRINT_PERSISTENCE_PAIRS std::cout << " [0," << get_diameter(e) << ")" << std::endl; #endif @@ -940,24 +932,23 @@ int main(int argc, char** argv) { columns_to_reduce.push_back(e); } std::reverse(columns_to_reduce.begin(), columns_to_reduce.end()); - + #ifdef PRINT_PERSISTENCE_PAIRS for (index_t i = 0; i < n; ++i) - if (dset.find(i) == i) - std::cout << " [0, )" << std::endl << std::flush; + if (dset.find(i) == i) std::cout << " [0, )" << std::endl << std::flush; #endif } - + for (index_t dim = 1; dim <= dim_max; ++dim) { rips_filtration_comparator comp(dist, dim + 1, binomial_coeff); rips_filtration_comparator comp_prev(dist, dim, binomial_coeff); - + hash_map pivot_column_index; pivot_column_index.reserve(columns_to_reduce.size()); - + compute_pairs(columns_to_reduce, pivot_column_index, dist, comp, comp_prev, dim, n, threshold, modulus, - multiplicative_inverse, binomial_coeff); - + multiplicative_inverse, binomial_coeff); + if (dim < dim_max) { assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, threshold, binomial_coeff); // std::cout << columns_to_reduce << std::endl; -- cgit v1.2.3 From 099ce879f4c729ec1fcbb148477aae64750bd90b Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 25 Aug 2016 23:17:00 +0200 Subject: packed entry_t --- ripser.cpp | 10 +++++++--- ripser.xcodeproj/project.pbxproj | 3 ++- 2 files changed, 9 insertions(+), 4 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index b4b528c..d3b969a 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -47,7 +47,7 @@ typedef float value_t; // typedef uint16_t value_t; typedef long index_t; -typedef long coefficient_t; +typedef short coefficient_t; class binomial_coeff_table { std::vector> B; @@ -132,18 +132,22 @@ std::vector vertices_of_simplex(const index_t simplex_index, const inde #ifdef USE_COEFFICIENTS struct entry_t { - index_t index; + index_t index : 8 * (sizeof(index_t) - sizeof(coefficient_t)); coefficient_t coefficient; entry_t(index_t _index, coefficient_t _coefficient) : index(_index), coefficient(_coefficient) {} entry_t(index_t _index) : index(_index), coefficient(1) {} entry_t() : index(0), coefficient(1) {} -}; +} __attribute__((packed)); + +static_assert(sizeof(entry_t) == sizeof(index_t), "size of entry_t is not the same as index_t"); + entry_t make_entry(index_t _index, coefficient_t _coefficient) { return entry_t(_index, _coefficient); } index_t get_index(entry_t e) { return e.index; } index_t get_coefficient(entry_t e) { return e.coefficient; } void set_coefficient(entry_t& e, const coefficient_t c) { e.coefficient = c; } + bool operator==(const entry_t& e1, const entry_t& e2) { return get_index(e1) == get_index(e2) && get_coefficient(e1) == get_coefficient(e2); } diff --git a/ripser.xcodeproj/project.pbxproj b/ripser.xcodeproj/project.pbxproj index 980df17..04bc3a6 100644 --- a/ripser.xcodeproj/project.pbxproj +++ b/ripser.xcodeproj/project.pbxproj @@ -209,9 +209,10 @@ buildSettings = { GCC_PREPROCESSOR_DEFINITIONS = ( "$(inherited)", - FILE_FORMAT_DISTANCE_MATRIX, + _FILE_FORMAT_DISTANCE_MATRIX, _FILE_FORMAT_LOWER_DISTANCE_MATRIX, _FILE_FORMAT_POINT_CLOUD, + USE_COEFFICIENTS, ); PRODUCT_NAME = "$(TARGET_NAME)"; }; -- cgit v1.2.3 From 76e9252dd33b97a438af8e05e6d4f0ffe0a625d0 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Sun, 28 Aug 2016 09:58:16 +0200 Subject: cleaned up detection of apparent persistence pairs --- ripser.cpp | 47 ++++++++++++++++++++--------------------------- 1 file changed, 20 insertions(+), 27 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index d3b969a..aec4e86 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -141,13 +141,11 @@ struct entry_t { static_assert(sizeof(entry_t) == sizeof(index_t), "size of entry_t is not the same as index_t"); - entry_t make_entry(index_t _index, coefficient_t _coefficient) { return entry_t(_index, _coefficient); } index_t get_index(entry_t e) { return e.index; } index_t get_coefficient(entry_t e) { return e.coefficient; } void set_coefficient(entry_t& e, const coefficient_t c) { e.coefficient = c; } - bool operator==(const entry_t& e1, const entry_t& e2) { return get_index(e1) == get_index(e2) && get_coefficient(e1) == get_coefficient(e2); } @@ -435,12 +433,12 @@ public: return entries.cbegin() + bounds[index]; } - template void append(Iterator begin, Iterator end) { + template void append_column(Iterator begin, Iterator end) { for (Iterator it = begin; it != end; ++it) { entries.push_back(*it); } bounds.push_back(entries.size()); } - void append() { bounds.push_back(entries.size()); } + void append_column() { bounds.push_back(entries.size()); } void push_back(ValueType e) { assert(0 < size()); @@ -454,8 +452,8 @@ public: --bounds.back(); } - template void append(const Collection collection) { - append(collection.cbegin(), collection.cend()); + template void append_column(const Collection collection) { + append_column(collection.cbegin(), collection.cend()); } }; @@ -527,7 +525,8 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map, - greater_diameter_or_smaller_index> working_coboundary; + greater_diameter_or_smaller_index> + working_coboundary; value_t diameter = get_diameter(column_to_reduce); @@ -539,18 +538,14 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map& columns_to_reduce, hash_map& columns_to_reduce, hash_map& columns_to_reduce, hash_map= 0); diameter_entry_t coface_entry = make_diameter_entry(coface_diameter, coface_index, coface_coefficient); coface_entries.push_back(coface_entry); - if (might_be_apparent_pair && (simplex_diameter == coface_diameter)) { + if (might_be_apparent_pair && (get_diameter(simplex) == coface_diameter)) { if (pivot_column_index.find(coface_index) == pivot_column_index.end()) { pivot = coface_entry; goto found_persistence_pair; @@ -629,7 +621,6 @@ void compute_pairs(std::vector& columns_to_reduce, hash_mapsecond; - column_to_add = columns_to_reduce[j]; continue; } } else { @@ -848,7 +839,9 @@ int main(int argc, char** argv) { for (index_t i = 1; i < argc; ++i) { const std::string arg(argv[i]); - if (arg == "--help") { print_usage_and_exit(0); } else if (arg == "--dim") { + if (arg == "--help") { + print_usage_and_exit(0); + } else if (arg == "--dim") { std::string parameter = std::string(argv[++i]); size_t next_pos; dim_max = std::stol(parameter, &next_pos); -- cgit v1.2.3 From 75f790a36a9dd6af3d1a6e39881f70474a55ece3 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Sun, 28 Aug 2016 09:23:19 +0200 Subject: added debug version to makefile --- Makefile | 7 +++++-- 1 file changed, 5 insertions(+), 2 deletions(-) diff --git a/Makefile b/Makefile index 331fdd6..73e09fb 100644 --- a/Makefile +++ b/Makefile @@ -1,7 +1,7 @@ build: ripser -all: ripser ripser-coeff ripser-reduction +all: ripser ripser-coeff ripser-reduction ripser-debug ripser: ripser.cpp @@ -13,6 +13,9 @@ ripser-coeff: ripser.cpp ripser-reduction: ripser.cpp c++ -std=c++11 ripser.cpp -o ripser-reduction -Ofast -D NDEBUG -D ASSEMBLE_REDUCTION_MATRIX +ripser-debug: ripser.cpp + c++ -std=c++11 ripser.cpp -o ripser-debug -g + clean: - rm -f ripser ripser-coeff ripser-reduction + rm -f ripser ripser-coeff ripser-reduction ripser-debug -- cgit v1.2.3 From 59c278afe11f99008b313ac48a7bf6bb5682c098 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Sun, 28 Aug 2016 10:03:52 +0200 Subject: changed to more general default format --- ripser.cpp | 2 +- 1 file changed, 1 insertion(+), 1 deletion(-) diff --git a/ripser.cpp b/ripser.cpp index aec4e86..5b4f053 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -826,7 +826,7 @@ int main(int argc, char** argv) { const char* filename = nullptr; - file_format format = LOWER_DISTANCE_MATRIX; + file_format format = DISTANCE_MATRIX; index_t dim_max = 1; value_t threshold = std::numeric_limits::max(); -- cgit v1.2.3 From 3f9fc342fef7a898b9242d39bd14ff4122c4036b Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Tue, 30 Aug 2016 18:19:58 +0200 Subject: cleanup --- ripser.cpp | 1 - 1 file changed, 1 deletion(-) diff --git a/ripser.cpp b/ripser.cpp index 5b4f053..5568f32 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -948,7 +948,6 @@ int main(int argc, char** argv) { if (dim < dim_max) { assemble_columns_to_reduce(columns_to_reduce, pivot_column_index, comp, dim, n, threshold, binomial_coeff); - // std::cout << columns_to_reduce << std::endl; } } } \ No newline at end of file -- cgit v1.2.3 From 3a1ff214fa33cba713f5a00ee13bfd0cebb0865f Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 15 Sep 2016 20:41:12 +0200 Subject: removed unneccesary check in file reader --- ripser.cpp | 2 +- 1 file changed, 1 insertion(+), 1 deletion(-) diff --git a/ripser.cpp b/ripser.cpp index 5568f32..9d8aa04 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -715,7 +715,7 @@ compressed_lower_distance_matrix read_point_cloud(std::istream& input_stream) { for (int i = 0; i < n; ++i) for (int j = 0; j < i; ++j) - if (i > j) distances.push_back(eucl_dist(i, j)); + distances.push_back(eucl_dist(i, j)); return compressed_lower_distance_matrix(std::move(distances)); } -- cgit v1.2.3 From 1ae58195ed01700ebe1cb0ece116be3c27a59faa Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 15 Sep 2016 20:45:16 +0200 Subject: pretty-print --- ripser.cpp | 7 +++---- 1 file changed, 3 insertions(+), 4 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index 9d8aa04..6241eb3 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -573,10 +573,10 @@ void compute_pairs(std::vector& columns_to_reduce, hash_map distances; for (int i = 0; i < n; ++i) - for (int j = 0; j < i; ++j) - distances.push_back(eucl_dist(i, j)); + for (int j = 0; j < i; ++j) distances.push_back(eucl_dist(i, j)); return compressed_lower_distance_matrix(std::move(distances)); } -- cgit v1.2.3 From 9f81b0f020b858fe95a10d35e984e37ff6d70175 Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 15 Sep 2016 20:53:13 +0200 Subject: explicit word sizes for index_t and coefficient_t --- ripser.cpp | 4 ++-- 1 file changed, 2 insertions(+), 2 deletions(-) diff --git a/ripser.cpp b/ripser.cpp index 6241eb3..600315b 100644 --- a/ripser.cpp +++ b/ripser.cpp @@ -46,8 +46,8 @@ template class hash_map : public std::unordered_map typedef float value_t; // typedef uint16_t value_t; -typedef long index_t; -typedef short coefficient_t; +typedef int64_t index_t; +typedef int16_t coefficient_t; class binomial_coeff_table { std::vector> B; -- cgit v1.2.3 From 7d5c7bb72b11f145190fa329949f478f987ceb4f Mon Sep 17 00:00:00 2001 From: Ulrich Bauer Date: Thu, 15 Sep 2016 21:09:07 +0200 Subject: updated version number in readme --- README.md | 2 +- 1 file changed, 1 insertion(+), 1 deletion(-) diff --git a/README.md b/README.md index 1f41a15..0a3f304 100644 --- a/README.md +++ b/README.md @@ -34,7 +34,7 @@ Ripser's efficiency is based on a few important concepts and principles: ### Version -[Latest release][latest-release]: 1.0 (July 2016) +[Latest release][latest-release]: 1.0.1 (September 2016) ### Building -- cgit v1.2.3