summaryrefslogtreecommitdiff
path: root/src/Persistence_representations/include
diff options
context:
space:
mode:
authormcarrier <mcarrier@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-02-16 15:43:29 +0000
committermcarrier <mcarrier@636b058d-ea47-450e-bf9e-a15bfbe3eedb>2018-02-16 15:43:29 +0000
commitff0dc023588e3b33bc4bc7f26ce1f68c647ae441 (patch)
treea6f839885acbbefe07ffeeca996eea77dc136e96 /src/Persistence_representations/include
parent69c683e663329d8410ca77c371f877bcc3bef906 (diff)
git-svn-id: svn+ssh://scm.gforge.inria.fr/svnroot/gudhi/branches/kernels@3251 636b058d-ea47-450e-bf9e-a15bfbe3eedb
Former-commit-id: 80f084fc990df6e5c6b60ac83514220aba2ceb5c
Diffstat (limited to 'src/Persistence_representations/include')
-rw-r--r--src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h78
-rw-r--r--src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h26
2 files changed, 59 insertions, 45 deletions
diff --git a/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h b/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h
index 2884885c..2b25b9a8 100644
--- a/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h
+++ b/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h
@@ -45,6 +45,7 @@
double pi = boost::math::constants::pi<double>();
using PD = std::vector<std::pair<double,double> >;
+using Weight = std::function<double (std::pair<double,double>) >;
namespace Gudhi {
namespace Persistence_representations {
@@ -53,11 +54,18 @@ class Persistence_weighted_gaussian{
protected:
PD diagram;
+ Weight weight;
+ double sigma;
+ int approx;
public:
- Persistence_weighted_gaussian(PD _diagram){diagram = _diagram;}
+ Persistence_weighted_gaussian(PD _diagram){diagram = _diagram; sigma = 1.0; approx = 1000; weight = arctan_weight;}
+ Persistence_weighted_gaussian(PD _diagram, double _sigma, int _approx, Weight _weight){diagram = _diagram; sigma = _sigma; approx = _approx; weight = _weight;}
PD get_diagram(){return this->diagram;}
+ double get_sigma(){return this->sigma;}
+ int get_approx(){return this->approx;}
+ Weight get_weight(){return this->weight;}
// **********************************
@@ -65,38 +73,37 @@ class Persistence_weighted_gaussian{
// **********************************
- static double pss_weight(std::pair<double,double> P){
- if(P.second > P.first) return 1;
+ static double pss_weight(std::pair<double,double> p){
+ if(p.second > p.first) return 1;
else return -1;
}
- static double arctan_weight(std::pair<double,double> P){
- return atan(P.second - P.first);
+ static double arctan_weight(std::pair<double,double> p){
+ return atan(p.second - p.first);
}
- template<class Weight = std::function<double (std::pair<double,double>) > >
- std::vector<std::pair<double,double> > Fourier_feat(PD D, std::vector<std::pair<double,double> > Z, Weight weight = arctan_weight){
- int m = D.size(); std::vector<std::pair<double,double> > B; int M = Z.size();
- for(int i = 0; i < M; i++){
- double d1 = 0; double d2 = 0; double zx = Z[i].first; double zy = Z[i].second;
- for(int j = 0; j < m; j++){
- double x = D[j].first; double y = D[j].second;
- d1 += weight(D[j])*cos(x*zx + y*zy);
- d2 += weight(D[j])*sin(x*zx + y*zy);
+ std::vector<std::pair<double,double> > Fourier_feat(PD diag, std::vector<std::pair<double,double> > z, Weight weight = arctan_weight){
+ int md = diag.size(); std::vector<std::pair<double,double> > b; int mz = z.size();
+ for(int i = 0; i < mz; i++){
+ double d1 = 0; double d2 = 0; double zx = z[i].first; double zy = z[i].second;
+ for(int j = 0; j < md; j++){
+ double x = diag[j].first; double y = diag[j].second;
+ d1 += weight(diag[j])*cos(x*zx + y*zy);
+ d2 += weight(diag[j])*sin(x*zx + y*zy);
}
- B.emplace_back(d1,d2);
+ b.emplace_back(d1,d2);
}
- return B;
+ return b;
}
- std::vector<std::pair<double,double> > random_Fourier(double sigma, int M = 1000){
- std::normal_distribution<double> distrib(0,1); std::vector<std::pair<double,double> > Z; std::random_device rd;
- for(int i = 0; i < M; i++){
+ std::vector<std::pair<double,double> > random_Fourier(double sigma, int m = 1000){
+ std::normal_distribution<double> distrib(0,1); std::vector<std::pair<double,double> > z; std::random_device rd;
+ for(int i = 0; i < m; i++){
std::mt19937 e1(rd()); std::mt19937 e2(rd());
double zx = distrib(e1); double zy = distrib(e2);
- Z.emplace_back(zx/sigma,zy/sigma);
+ z.emplace_back(zx/sigma,zy/sigma);
}
- return Z;
+ return z;
}
@@ -106,32 +113,33 @@ class Persistence_weighted_gaussian{
// **********************************
- template<class Weight = std::function<double (std::pair<double,double>) > >
- double compute_scalar_product(Persistence_weighted_gaussian second, double sigma, Weight weight = arctan_weight, int m = 1000){
+ double compute_scalar_product(Persistence_weighted_gaussian second){
PD diagram1 = this->diagram; PD diagram2 = second.diagram;
- if(m == -1){
+ if(this->approx == -1){
int num_pts1 = diagram1.size(); int num_pts2 = diagram2.size(); double k = 0;
for(int i = 0; i < num_pts1; i++)
for(int j = 0; j < num_pts2; j++)
- k += weight(diagram1[i])*weight(diagram2[j])*exp(-((diagram1[i].first - diagram2[j].first) * (diagram1[i].first - diagram2[j].first) +
- (diagram1[i].second - diagram2[j].second) * (diagram1[i].second - diagram2[j].second))
- /(2*sigma*sigma));
+ k += this->weight(diagram1[i])*this->weight(diagram2[j])*exp(-((diagram1[i].first - diagram2[j].first) * (diagram1[i].first - diagram2[j].first) +
+ (diagram1[i].second - diagram2[j].second) * (diagram1[i].second - diagram2[j].second))
+ /(2*this->sigma*this->sigma));
return k;
}
else{
- std::vector<std::pair<double,double> > z = random_Fourier(sigma, m);
- std::vector<std::pair<double,double> > b1 = Fourier_feat(diagram1,z,weight);
- std::vector<std::pair<double,double> > b2 = Fourier_feat(diagram2,z,weight);
- double d = 0; for(int i = 0; i < m; i++) d += b1[i].first*b2[i].first + b1[i].second*b2[i].second;
- return d/m;
+ std::vector<std::pair<double,double> > z = random_Fourier(this->sigma, this->approx);
+ std::vector<std::pair<double,double> > b1 = Fourier_feat(diagram1,z,this->weight);
+ std::vector<std::pair<double,double> > b2 = Fourier_feat(diagram2,z,this->weight);
+ double d = 0; for(int i = 0; i < this->approx; i++) d += b1[i].first*b2[i].first + b1[i].second*b2[i].second;
+ return d/this->approx;
}
}
- template<class Weight = std::function<double (std::pair<double,double>) > >
- double distance(Persistence_weighted_gaussian second, double sigma, Weight weight = arctan_weight, int m = 1000, double power = 1) {
- return std::pow(this->compute_scalar_product(*this, sigma, weight, m) + second.compute_scalar_product(second, sigma, weight, m)-2*this->compute_scalar_product(second, sigma, weight, m), power/2.0);
+ double distance(Persistence_weighted_gaussian second, double power = 1) {
+ if(this->sigma != second.get_sigma() || this->approx != second.get_approx()){
+ std::cout << "Error: different representations!" << std::endl; return 0;
+ }
+ else return std::pow(this->compute_scalar_product(*this) + second.compute_scalar_product(second)-2*this->compute_scalar_product(second), power/2.0);
}
diff --git a/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h b/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h
index 4fa6151f..ad1a6c42 100644
--- a/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h
+++ b/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h
@@ -53,11 +53,16 @@ class Sliced_Wasserstein {
protected:
PD diagram;
+ int approx;
+ double sigma;
public:
- Sliced_Wasserstein(PD _diagram){diagram = _diagram;}
+ Sliced_Wasserstein(PD _diagram){diagram = _diagram; approx = 100; sigma = 0.001;}
+ Sliced_Wasserstein(PD _diagram, double _sigma, int _approx){diagram = _diagram; approx = _approx; sigma = _sigma;}
PD get_diagram(){return this->diagram;}
+ int get_approx(){return this->approx;}
+ double get_sigma(){return this->sigma;}
// **********************************
@@ -130,11 +135,11 @@ class Sliced_Wasserstein {
// Scalar product + distance.
// **********************************
- double compute_sliced_wasserstein_distance(Sliced_Wasserstein second, int approx) {
+ double compute_sliced_wasserstein_distance(Sliced_Wasserstein second) {
PD diagram1 = this->diagram; PD diagram2 = second.diagram; double sw = 0;
- if(approx == -1){
+ if(this->approx == -1){
// Add projections onto diagonal.
int n1, n2; n1 = diagram1.size(); n2 = diagram2.size(); double max_ordinate = std::numeric_limits<double>::lowest();
@@ -226,7 +231,7 @@ class Sliced_Wasserstein {
else{
- double step = pi/approx;
+ double step = pi/this->approx;
// Add projections onto diagonal.
int n1, n2; n1 = diagram1.size(); n2 = diagram2.size();
@@ -238,7 +243,7 @@ class Sliced_Wasserstein {
// Sort and compare all projections.
#ifdef GUDHI_USE_TBB
- tbb::parallel_for(0, approx, [&](int i){
+ tbb::parallel_for(0, this->approx, [&](int i){
std::vector<std::pair<int,double> > l1, l2;
for (int j = 0; j < n; j++){
l1.emplace_back( j, diagram1[j].first*cos(-pi/2+i*step) + diagram1[j].second*sin(-pi/2+i*step) );
@@ -250,7 +255,7 @@ class Sliced_Wasserstein {
sw += f*step;
});
#else
- for (int i = 0; i < approx; i++){
+ for (int i = 0; i < this->approx; i++){
std::vector<std::pair<int,double> > l1, l2;
for (int j = 0; j < n; j++){
l1.emplace_back( j, diagram1[j].first*cos(-pi/2+i*step) + diagram1[j].second*sin(-pi/2+i*step) );
@@ -268,12 +273,13 @@ class Sliced_Wasserstein {
}
- double compute_scalar_product(Sliced_Wasserstein second, double sigma, int approx = 100) {
- return std::exp(-compute_sliced_wasserstein_distance(second, approx)/(2*sigma*sigma));
+ double compute_scalar_product(Sliced_Wasserstein second){
+ return std::exp(-compute_sliced_wasserstein_distance(second)/(2*this->sigma*this->sigma));
}
- double distance(Sliced_Wasserstein second, double sigma, int approx = 100, double power = 1) {
- return std::pow(this->compute_scalar_product(*this, sigma, approx) + second.compute_scalar_product(second, sigma, approx)-2*this->compute_scalar_product(second, sigma, approx), power/2.0);
+ double distance(Sliced_Wasserstein second, double power = 1) {
+ if(this->sigma != second.sigma || this->approx != second.approx){std::cout << "Error: different representations!" << std::endl; return 0;}
+ else return std::pow(this->compute_scalar_product(*this) + second.compute_scalar_product(second)-2*this->compute_scalar_product(second), power/2.0);
}