summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorUlrich Bauer <mail@ulrich-bauer.org>2019-06-24 19:24:34 +0200
committerUlrich Bauer <mail@ulrich-bauer.org>2019-06-24 19:24:34 +0200
commitcb6163df1dd96e2f65c5b4506fd88c1266562e86 (patch)
tree69f8bfa7a7d9a5629cd6d749a70920e7c9de5199
parent246f1fe3efcb794b88c95865477be2b7dc0cf709 (diff)
structs instead of classes
-rw-r--r--README.md21
-rw-r--r--ripser.cpp62
2 files changed, 33 insertions, 50 deletions
diff --git a/README.md b/README.md
index c8f82cc..3c556c0 100644
--- a/README.md
+++ b/README.md
@@ -1,6 +1,6 @@
# Ripser
-Copyright © 2015–2018 [Ulrich Bauer].
+Copyright © 2015–2019 [Ulrich Bauer].
### Description
@@ -24,6 +24,8 @@ Input formats currently supported by Ripser:
- comma-separated values upper triangular distance matrix (MATLAB output from the function `pdist`)
- comma-separated values full distance matrix
- [DIPHA] distance matrix data
+ - sparse distance matrix in Sparse Triplet format
+ - binary lower triangular distance matrix
- point cloud data
Ripser's efficiency is based on a few important concepts and principles, building on key previous and concurrent developments by other researchers in computational topology:
@@ -37,7 +39,7 @@ Ripser's efficiency is based on a few important concepts and principles, buildin
### Version
-[Latest release][latest-release]: 1.0.1 (September 2016)
+[Latest release][latest-release]: 1.1 (June 2019)
### Building
@@ -56,16 +58,14 @@ 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 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)
- `USE_GOOGLE_HASHMAP`: enable support for Google's [sparsehash] data structure; may further reduce memory footprint
-For example, to build Ripser with support for coefficients:
+For example, to build Ripser with support for Google's hashmap:
```sh
-$ c++ -std=c++11 ripser.cpp -o ripser -Ofast -D NDEBUG -D USE_COEFFICIENTS
+$ c++ -std=c++11 ripser.cpp -o ripser -Ofast -D NDEBUG -D USE_GOOGLE_HASHMAP
```
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.
@@ -78,19 +78,20 @@ The input is given either in a file whose name is passed as an argument, or thro
- `distance`: full distance matrix; similar to the above, but for all entries of the distance matrix. One line per row of the matrix; only the part below the diagonal is actually read.
- `dipha`: DIPHA distance matrix as described on the [DIPHA] website.
- `point-cloud`: point cloud; a comma (or whitespace, or other non-numerical character) separated list of coordinates of the points in some Euclidean space, one point per line.
+ - `sparse`: sparse distance matrix in Sparse Triplet format
+ - `binary`: lower distance matrix in binary file format; a sequence of the distance matrix entries below the diagonal in 64 bit double format (IEEE 754, little endian).
- `--dim k`: compute persistent homology up to dimension *k*.
- `--threshold t`: compute Rips complexes up to diameter *t*.
- - `--modulus p`: compute homology with coefficients in the prime field Z/*p*Z (only available when built with the option `USE_COEFFICIENTS`).
+ - `--modulus p`: compute homology with coefficients in the prime field Z/*p*Z.
+ - `--ratio <r>`: only show persistence pairs with death/birth ratio > *r*.
-
-### Planned features
+### Experimental features
The following features are currently planned for future versions:
- computation of representative cycles for persistent homology (currenly only *co*cycles are computed)
- - support for sparse distance matrices
Prototype implementations are already avaliable; please contact the author if one of these features might be relevant for your research.
diff --git a/ripser.cpp b/ripser.cpp
index 3d08300..bdb452a 100644
--- a/ripser.cpp
+++ b/ripser.cpp
@@ -41,9 +41,6 @@
//#define USE_GOOGLE_HASHMAP
-#include <algorithm>
-#include <cassert>
-#include <cmath>
#include <fstream>
#include <iostream>
#include <numeric>
@@ -141,8 +138,7 @@ typedef std::pair<index_t, value_t> index_diameter_t;
index_t get_index(const index_diameter_t& i) { return i.first; }
value_t get_diameter(const index_diameter_t& i) { return i.second; }
-class diameter_entry_t : public std::pair<value_t, entry_t> {
-public:
+struct diameter_entry_t : std::pair<value_t, entry_t> {
diameter_entry_t() {}
diameter_entry_t(value_t _diameter, index_t _index, coefficient_t _coefficient)
: std::pair<value_t, entry_t>(_diameter, make_entry(_index, _coefficient)) {}
@@ -172,8 +168,7 @@ template <typename Entry> struct greater_diameter_or_smaller_index {
enum compressed_matrix_layout { LOWER_TRIANGULAR, UPPER_TRIANGULAR };
-template <compressed_matrix_layout Layout> class compressed_distance_matrix {
-public:
+template <compressed_matrix_layout Layout> struct compressed_distance_matrix {
std::vector<value_t> distances;
std::vector<value_t*> rows;
@@ -199,8 +194,7 @@ public:
size_t size() const { return rows.size(); }
};
-class sparse_distance_matrix {
-public:
+struct sparse_distance_matrix {
std::vector<std::vector<index_diameter_t>> neighbors;
index_t num_edges;
@@ -258,8 +252,7 @@ value_t compressed_distance_matrix<LOWER_TRIANGULAR>::operator()(const index_t i
typedef compressed_distance_matrix<LOWER_TRIANGULAR> compressed_lower_distance_matrix;
typedef compressed_distance_matrix<UPPER_TRIANGULAR> compressed_upper_distance_matrix;
-class euclidean_distance_matrix {
-public:
+struct euclidean_distance_matrix {
std::vector<std::vector<value_t>> points;
euclidean_distance_matrix(std::vector<std::vector<value_t>>&& _points)
@@ -393,7 +386,6 @@ template <typename DistanceMatrix> class ripser {
coefficient_t modulus;
const binomial_coeff_table binomial_coeff;
std::vector<coefficient_t> multiplicative_inverse;
- mutable std::vector<index_t> vertices;
mutable std::vector<diameter_entry_t> coface_entries;
public:
@@ -474,7 +466,6 @@ public:
union_find dset(n);
edges = get_edges();
-
std::sort(edges.rbegin(), edges.rend(),
greater_diameter_or_smaller_index<diameter_index_t>());
@@ -553,11 +544,7 @@ public:
std::priority_queue<diameter_entry_t, std::vector<diameter_entry_t>,
greater_diameter_or_smaller_index<diameter_entry_t>>
- working_reduction_column;
-
- std::priority_queue<diameter_entry_t, std::vector<diameter_entry_t>,
- greater_diameter_or_smaller_index<diameter_entry_t>>
- working_coboundary;
+ working_reduction_column, working_coboundary;
value_t diameter = get_diameter(column_to_reduce);
@@ -644,7 +631,6 @@ public:
std::vector<diameter_index_t> get_edges();
void compute_barcodes() {
-
std::vector<diameter_index_t> simplices, columns_to_reduce;
compute_dim_0_pairs(simplices, columns_to_reduce);
@@ -655,16 +641,14 @@ public:
compute_pairs(columns_to_reduce, pivot_column_index, dim);
- if (dim < dim_max) {
+ if (dim < dim_max)
assemble_columns_to_reduce(simplices, columns_to_reduce, pivot_column_index,
dim + 1);
- }
}
}
};
template <> class ripser<compressed_lower_distance_matrix>::simplex_coboundary_enumerator {
-private:
index_t idx_below, idx_above, v, k;
std::vector<index_t> vertices;
const diameter_entry_t simplex;
@@ -704,16 +688,15 @@ public:
};
template <> class ripser<sparse_distance_matrix>::simplex_coboundary_enumerator {
-private:
const ripser& parent;
index_t idx_below, idx_above, v, k, max_vertex_below;
+ std::vector<index_t> vertices;
const diameter_entry_t simplex;
const coefficient_t modulus;
const sparse_distance_matrix& dist;
const binomial_coeff_table& binomial_coeff;
- std::vector<index_t> vertices;
std::vector<std::vector<index_diameter_t>::const_reverse_iterator>& neighbor_it;
std::vector<std::vector<index_diameter_t>::const_reverse_iterator>& neighbor_end;
index_diameter_t x;
@@ -800,7 +783,7 @@ enum file_format {
POINT_CLOUD,
DIPHA,
SPARSE,
- RIPSER
+ BINARY
};
template <typename T> T read(std::istream& s) {
@@ -841,9 +824,7 @@ compressed_lower_distance_matrix read_point_cloud(std::istream& input_stream) {
}
sparse_distance_matrix read_sparse_distance_matrix(std::istream& input_stream) {
-
std::vector<std::vector<index_diameter_t>> neighbors;
-
index_t num_edges = 0;
std::string line;
@@ -931,7 +912,7 @@ compressed_lower_distance_matrix read_dipha(std::istream& input_stream) {
return compressed_lower_distance_matrix(std::move(distances));
}
-compressed_lower_distance_matrix read_ripser(std::istream& input_stream) {
+compressed_lower_distance_matrix read_binary(std::istream& input_stream) {
std::vector<value_t> distances;
while (!input_stream.eof()) distances.push_back(read<value_t>(input_stream));
return compressed_lower_distance_matrix(std::move(distances));
@@ -950,7 +931,7 @@ compressed_lower_distance_matrix read_file(std::istream& input_stream, file_form
case DIPHA:
return read_dipha(input_stream);
default:
- return read_ripser(input_stream);
+ return read_binary(input_stream);
}
}
@@ -973,18 +954,19 @@ void print_usage_and_exit(int exit_code) {
<< " dipha (distance matrix in DIPHA file format)" << std::endl
<< " sparse (sparse distance matrix in Sparse Triplet format)"
<< std::endl
- << " ripser (distance matrix in Ripser binary file format)"
+ << " binary (lower triangular distance matrix in binary format)"
<< std::endl
- << " --dim <k> compute persistent homology up to dimension <k>" << std::endl
- << " --threshold <t> compute Rips complexes up to diameter <t>" << std::endl
- << " --modulus <p> compute homology with coefficients in the prime field Z/<p>Z"
+ << " --dim <k> compute persistent homology up to dimension k" << std::endl
+ << " --threshold <t> compute Rips complexes up to diameter t" << std::endl
+ << " --modulus <p> compute homology with coefficients in the prime field Z/pZ"
+ << std::endl
+ << " --ratio <r> only show persistence pairs with death/birth ratio > r"
<< std::endl;
exit(exit_code);
}
int main(int argc, char** argv) {
-
const char* filename = nullptr;
file_format format = DISTANCE_MATRIX;
@@ -1017,20 +999,20 @@ int main(int argc, char** argv) {
if (next_pos != parameter.size()) print_usage_and_exit(-1);
} else if (arg == "--format") {
std::string parameter = std::string(argv[++i]);
- if (parameter == "lower-distance")
+ if (parameter.rfind("lower", 0) == 0)
format = LOWER_DISTANCE_MATRIX;
- else if (parameter == "upper-distance")
+ else if (parameter.rfind("upper", 0) == 0)
format = UPPER_DISTANCE_MATRIX;
- else if (parameter == "distance")
+ else if (parameter.rfind("dist", 0) == 0)
format = DISTANCE_MATRIX;
- else if (parameter == "point-cloud")
+ else if (parameter.rfind("point", 0) == 0)
format = POINT_CLOUD;
else if (parameter == "dipha")
format = DIPHA;
else if (parameter == "sparse")
format = SPARSE;
- else if (parameter == "ripser")
- format = RIPSER;
+ else if (parameter == "binary")
+ format = BINARY;
else
print_usage_and_exit(-1);
} else if (arg == "--modulus") {