summaryrefslogtreecommitdiff
path: root/src/Persistence_representations
diff options
context:
space:
mode:
Diffstat (limited to 'src/Persistence_representations')
-rw-r--r--src/Persistence_representations/example/persistence_weighted_gaussian.cpp96
-rw-r--r--src/Persistence_representations/example/sliced_wasserstein.cpp55
-rw-r--r--src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h143
-rw-r--r--src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h285
4 files changed, 579 insertions, 0 deletions
diff --git a/src/Persistence_representations/example/persistence_weighted_gaussian.cpp b/src/Persistence_representations/example/persistence_weighted_gaussian.cpp
new file mode 100644
index 00000000..e95b9445
--- /dev/null
+++ b/src/Persistence_representations/example/persistence_weighted_gaussian.cpp
@@ -0,0 +1,96 @@
+/* This file is part of the Gudhi Library. The Gudhi library
+ * (Geometric Understanding in Higher Dimensions) is a generic C++
+ * library for computational topology.
+ *
+ * Author(s): Mathieu Carriere
+ *
+ * Copyright (C) 2018 INRIA (France)
+ *
+ * This program is free software: you can redistribute it and/or modify
+ * it under the terms of the GNU General Public License as published by
+ * the Free Software Foundation, either version 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ * GNU General Public License for more details.
+ *
+ * You should have received a copy of the GNU General Public License
+ * along with this program. If not, see <http://www.gnu.org/licenses/>.
+ */
+
+#include <gudhi/Persistence_weighted_gaussian.h>
+
+#include <iostream>
+#include <vector>
+#include <utility>
+
+using PD = std::vector<std::pair<double,double> >;
+using PWG = Gudhi::Persistence_representations::Persistence_weighted_gaussian;
+
+int main(int argc, char** argv) {
+
+ std::vector<std::pair<double, double> > persistence1;
+ std::vector<std::pair<double, double> > persistence2;
+
+ persistence1.push_back(std::make_pair(1, 2));
+ persistence1.push_back(std::make_pair(6, 8));
+ persistence1.push_back(std::make_pair(0, 4));
+ persistence1.push_back(std::make_pair(3, 8));
+
+ persistence2.push_back(std::make_pair(2, 9));
+ persistence2.push_back(std::make_pair(1, 6));
+ persistence2.push_back(std::make_pair(3, 5));
+ persistence2.push_back(std::make_pair(6, 10));
+
+ PWG PWG1(persistence1);
+ PWG PWG2(persistence2);
+ double sigma = 1;
+ double tau = 1;
+ int m = 1000;
+
+
+
+ // Linear PWG
+
+ std::cout << PWG1.compute_scalar_product (PWG2, sigma, PWG::arctan_weight, m) << std::endl;
+ std::cout << PWG1.compute_scalar_product (PWG2, sigma, PWG::arctan_weight, -1) << std::endl;
+
+ std::cout << PWG1.distance (PWG2, sigma, PWG::arctan_weight, m) << std::endl;
+ std::cout << PWG1.distance (PWG2, sigma, PWG::arctan_weight, -1) << std::endl;
+
+
+
+
+
+
+
+ // Gaussian PWG
+
+ std::cout << std::exp( -PWG1.distance (PWG2, sigma, PWG::arctan_weight, m, 2) ) / (2*tau*tau) << std::endl;
+ std::cout << std::exp( -PWG1.distance (PWG2, sigma, PWG::arctan_weight, -1, 2) ) / (2*tau*tau) << std::endl;
+
+
+
+
+
+
+
+ // PSS
+
+ PD pd1 = persistence1; int numpts = persistence1.size(); for(int i = 0; i < numpts; i++) pd1.emplace_back(persistence1[i].second,persistence1[i].first);
+ PD pd2 = persistence2; numpts = persistence2.size(); for(int i = 0; i < numpts; i++) pd2.emplace_back(persistence2[i].second,persistence2[i].first);
+
+ PWG pwg1(pd1);
+ PWG pwg2(pd2);
+
+ std::cout << pwg1.compute_scalar_product (pwg2, 2*std::sqrt(sigma), PWG::pss_weight, m) / (16*pi*sigma) << std::endl;
+ std::cout << pwg1.compute_scalar_product (pwg2, 2*std::sqrt(sigma), PWG::pss_weight, -1) / (16*pi*sigma) << std::endl;
+
+ std::cout << pwg1.distance (pwg2, 2*std::sqrt(sigma), PWG::pss_weight, m) / (16*pi*sigma) << std::endl;
+ std::cout << pwg1.distance (pwg2, 2*std::sqrt(sigma), PWG::pss_weight, -1) / (16*pi*sigma) << std::endl;
+
+
+ return 0;
+}
diff --git a/src/Persistence_representations/example/sliced_wasserstein.cpp b/src/Persistence_representations/example/sliced_wasserstein.cpp
new file mode 100644
index 00000000..673d8474
--- /dev/null
+++ b/src/Persistence_representations/example/sliced_wasserstein.cpp
@@ -0,0 +1,55 @@
+/* This file is part of the Gudhi Library. The Gudhi library
+ * (Geometric Understanding in Higher Dimensions) is a generic C++
+ * library for computational topology.
+ *
+ * Author(s): Mathieu Carriere
+ *
+ * Copyright (C) 2018 INRIA (France)
+ *
+ * This program is free software: you can redistribute it and/or modify
+ * it under the terms of the GNU General Public License as published by
+ * the Free Software Foundation, either version 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ * GNU General Public License for more details.
+ *
+ * You should have received a copy of the GNU General Public License
+ * along with this program. If not, see <http://www.gnu.org/licenses/>.
+ */
+
+#include <gudhi/Sliced_Wasserstein.h>
+
+#include <iostream>
+#include <vector>
+#include <utility>
+
+using SW = Gudhi::Persistence_representations::Sliced_Wasserstein;
+
+int main(int argc, char** argv) {
+
+ std::vector<std::pair<double, double> > persistence1;
+ std::vector<std::pair<double, double> > persistence2;
+
+ persistence1.push_back(std::make_pair(1, 2));
+ persistence1.push_back(std::make_pair(6, 8));
+ persistence1.push_back(std::make_pair(0, 4));
+ persistence1.push_back(std::make_pair(3, 8));
+
+ persistence2.push_back(std::make_pair(2, 9));
+ persistence2.push_back(std::make_pair(1, 6));
+ persistence2.push_back(std::make_pair(3, 5));
+ persistence2.push_back(std::make_pair(6, 10));
+
+ SW SW1(persistence1);
+ SW SW2(persistence2);
+
+ std::cout << SW1.compute_sliced_wasserstein_distance(SW2,100) << std::endl;
+ std::cout << SW1.compute_sliced_wasserstein_distance(SW2,-1) << std::endl;
+ std::cout << SW1.compute_scalar_product(SW2,1,100) << std::endl;
+ std::cout << SW1.distance(SW2,1,100,1) << std::endl;
+
+ return 0;
+}
diff --git a/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h b/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h
new file mode 100644
index 00000000..2884885c
--- /dev/null
+++ b/src/Persistence_representations/include/gudhi/Persistence_weighted_gaussian.h
@@ -0,0 +1,143 @@
+/* This file is part of the Gudhi Library. The Gudhi library
+ * (Geometric Understanding in Higher Dimensions) is a generic C++
+ * library for computational topology.
+ *
+ * Author(s): Mathieu Carriere
+ *
+ * Copyright (C) 2018 INRIA (France)
+ *
+ * This program is free software: you can redistribute it and/or modify
+ * it under the terms of the GNU General Public License as published by
+ * the Free Software Foundation, either version 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ * GNU General Public License for more details.
+ *
+ * You should have received a copy of the GNU General Public License
+ * along with this program. If not, see <http://www.gnu.org/licenses/>.
+ */
+
+#ifndef PERSISTENCE_WEIGHTED_GAUSSIAN_H_
+#define PERSISTENCE_WEIGHTED_GAUSSIAN_H_
+
+#ifdef GUDHI_USE_TBB
+#include <tbb/parallel_for.h>
+#endif
+
+// gudhi include
+#include <gudhi/read_persistence_from_file.h>
+
+// standard include
+#include <cmath>
+#include <iostream>
+#include <vector>
+#include <limits>
+#include <fstream>
+#include <sstream>
+#include <algorithm>
+#include <string>
+#include <utility>
+#include <functional>
+#include <boost/math/constants/constants.hpp>
+
+double pi = boost::math::constants::pi<double>();
+using PD = std::vector<std::pair<double,double> >;
+
+namespace Gudhi {
+namespace Persistence_representations {
+
+class Persistence_weighted_gaussian{
+
+ protected:
+ PD diagram;
+
+ public:
+
+ Persistence_weighted_gaussian(PD _diagram){diagram = _diagram;}
+ PD get_diagram(){return this->diagram;}
+
+
+ // **********************************
+ // Utils.
+ // **********************************
+
+
+ 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);
+ }
+
+ 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);
+ }
+ B.emplace_back(d1,d2);
+ }
+ 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::mt19937 e1(rd()); std::mt19937 e2(rd());
+ double zx = distrib(e1); double zy = distrib(e2);
+ Z.emplace_back(zx/sigma,zy/sigma);
+ }
+ return Z;
+ }
+
+
+
+ // **********************************
+ // Scalar product + distance.
+ // **********************************
+
+
+ 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){
+
+ PD diagram1 = this->diagram; PD diagram2 = second.diagram;
+
+ if(m == -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));
+ 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;
+ }
+ }
+
+ 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);
+ }
+
+
+};
+
+} // namespace Persistence_weighted_gaussian
+} // namespace Gudhi
+
+#endif // PERSISTENCE_WEIGHTED_GAUSSIAN_H_
diff --git a/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h b/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h
new file mode 100644
index 00000000..4fa6151f
--- /dev/null
+++ b/src/Persistence_representations/include/gudhi/Sliced_Wasserstein.h
@@ -0,0 +1,285 @@
+/* This file is part of the Gudhi Library. The Gudhi library
+ * (Geometric Understanding in Higher Dimensions) is a generic C++
+ * library for computational topology.
+ *
+ * Author(s): Mathieu Carriere
+ *
+ * Copyright (C) 2018 INRIA (France)
+ *
+ * This program is free software: you can redistribute it and/or modify
+ * it under the terms of the GNU General Public License as published by
+ * the Free Software Foundation, either version 3 of the License, or
+ * (at your option) any later version.
+ *
+ * This program is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
+ * GNU General Public License for more details.
+ *
+ * You should have received a copy of the GNU General Public License
+ * along with this program. If not, see <http://www.gnu.org/licenses/>.
+ */
+
+#ifndef SLICED_WASSERSTEIN_H_
+#define SLICED_WASSERSTEIN_H_
+
+#ifdef GUDHI_USE_TBB
+#include <tbb/parallel_for.h>
+#endif
+
+// gudhi include
+#include <gudhi/read_persistence_from_file.h>
+
+// standard include
+#include <cmath>
+#include <iostream>
+#include <vector>
+#include <limits>
+#include <fstream>
+#include <sstream>
+#include <algorithm>
+#include <string>
+#include <utility>
+#include <functional>
+#include <boost/math/constants/constants.hpp>
+
+double pi = boost::math::constants::pi<double>();
+using PD = std::vector<std::pair<double,double> >;
+
+namespace Gudhi {
+namespace Persistence_representations {
+
+class Sliced_Wasserstein {
+
+ protected:
+ PD diagram;
+
+ public:
+
+ Sliced_Wasserstein(PD _diagram){diagram = _diagram;}
+ PD get_diagram(){return this->diagram;}
+
+
+ // **********************************
+ // Utils.
+ // **********************************
+
+ // Compute the angle formed by two points of a PD
+ double compute_angle(PD diag, int i, int j){
+ std::pair<double,double> vect; double x1,y1, x2,y2;
+ x1 = diag[i].first; y1 = diag[i].second;
+ x2 = diag[j].first; y2 = diag[j].second;
+ if (y1 - y2 > 0){
+ vect.first = y1 - y2;
+ vect.second = x2 - x1;}
+ else{
+ if(y1 - y2 < 0){
+ vect.first = y2 - y1;
+ vect.second = x1 - x2;
+ }
+ else{
+ vect.first = 0;
+ vect.second = abs(x1 - x2);}
+ }
+ double norm = std::sqrt(vect.first*vect.first + vect.second*vect.second);
+ return asin(vect.second/norm);
+ }
+
+ // Compute the integral of |cos()| between alpha and beta, valid only if alpha is in [-pi,pi] and beta-alpha is in [0,pi]
+ double compute_int_cos(const double & alpha, const double & beta){
+ double res = 0;
+ if (alpha >= 0 && alpha <= pi){
+ if (cos(alpha) >= 0){
+ if(pi/2 <= beta){res = 2-sin(alpha)-sin(beta);}
+ else{res = sin(beta)-sin(alpha);}
+ }
+ else{
+ if(1.5*pi <= beta){res = 2+sin(alpha)+sin(beta);}
+ else{res = sin(alpha)-sin(beta);}
+ }
+ }
+ if (alpha >= -pi && alpha <= 0){
+ if (cos(alpha) <= 0){
+ if(-pi/2 <= beta){res = 2+sin(alpha)+sin(beta);}
+ else{res = sin(alpha)-sin(beta);}
+ }
+ else{
+ if(pi/2 <= beta){res = 2-sin(alpha)-sin(beta);}
+ else{res = sin(beta)-sin(alpha);}
+ }
+ }
+ return res;
+ }
+
+ double compute_int(const double & theta1, const double & theta2, const int & p, const int & q, const PD & PD1, const PD & PD2){
+ double norm = std::sqrt( (PD1[p].first-PD2[q].first)*(PD1[p].first-PD2[q].first) + (PD1[p].second-PD2[q].second)*(PD1[p].second-PD2[q].second) );
+ double angle1;
+ if (PD1[p].first > PD2[q].first)
+ angle1 = theta1 - asin( (PD1[p].second-PD2[q].second)/norm );
+ else
+ angle1 = theta1 - asin( (PD2[q].second-PD1[p].second)/norm );
+ double angle2 = angle1 + theta2 - theta1;
+ double integral = compute_int_cos(angle1,angle2);
+ return norm*integral;
+ }
+
+
+
+
+ // **********************************
+ // Scalar product + distance.
+ // **********************************
+
+ double compute_sliced_wasserstein_distance(Sliced_Wasserstein second, int approx) {
+
+ PD diagram1 = this->diagram; PD diagram2 = second.diagram; double sw = 0;
+
+ if(approx == -1){
+
+ // Add projections onto diagonal.
+ int n1, n2; n1 = diagram1.size(); n2 = diagram2.size(); double max_ordinate = std::numeric_limits<double>::lowest();
+ for (int i = 0; i < n2; i++){
+ max_ordinate = std::max(max_ordinate, diagram2[i].second);
+ diagram1.emplace_back( (diagram2[i].first+diagram2[i].second)/2, (diagram2[i].first+diagram2[i].second)/2 );
+ }
+ for (int i = 0; i < n1; i++){
+ max_ordinate = std::max(max_ordinate, diagram1[i].second);
+ diagram2.emplace_back( (diagram1[i].first+diagram1[i].second)/2, (diagram1[i].first+diagram1[i].second)/2 );
+ }
+ int num_pts_dgm = diagram1.size();
+
+ // Slightly perturb the points so that the PDs are in generic positions.
+ int mag = 0; while(max_ordinate > 10){mag++; max_ordinate/=10;}
+ double thresh = pow(10,-5+mag);
+ srand(time(NULL));
+ for (int i = 0; i < num_pts_dgm; i++){
+ diagram1[i].first += thresh*(1.0-2.0*rand()/RAND_MAX); diagram1[i].second += thresh*(1.0-2.0*rand()/RAND_MAX);
+ diagram2[i].first += thresh*(1.0-2.0*rand()/RAND_MAX); diagram2[i].second += thresh*(1.0-2.0*rand()/RAND_MAX);
+ }
+
+ // Compute all angles in both PDs.
+ std::vector<std::pair<double, std::pair<int,int> > > angles1, angles2;
+ for (int i = 0; i < num_pts_dgm; i++){
+ for (int j = i+1; j < num_pts_dgm; j++){
+ double theta1 = compute_angle(diagram1,i,j); double theta2 = compute_angle(diagram2,i,j);
+ angles1.emplace_back(theta1, std::pair<int,int>(i,j));
+ angles2.emplace_back(theta2, std::pair<int,int>(i,j));
+ }
+ }
+
+ // Sort angles.
+ std::sort(angles1.begin(), angles1.end(), [=](std::pair<double, std::pair<int,int> >& p1, const std::pair<double, std::pair<int,int> >& p2){return (p1.first < p2.first);});
+ std::sort(angles2.begin(), angles2.end(), [=](std::pair<double, std::pair<int,int> >& p1, const std::pair<double, std::pair<int,int> >& p2){return (p1.first < p2.first);});
+
+ // Initialize orders of the points of both PDs (given by ordinates when theta = -pi/2).
+ std::vector<int> orderp1, orderp2;
+ for (int i = 0; i < num_pts_dgm; i++){ orderp1.push_back(i); orderp2.push_back(i); }
+ std::sort( orderp1.begin(), orderp1.end(), [=](int i, int j){ if(diagram1[i].second != diagram1[j].second) return (diagram1[i].second < diagram1[j].second); else return (diagram1[i].first > diagram1[j].first); } );
+ std::sort( orderp2.begin(), orderp2.end(), [=](int i, int j){ if(diagram2[i].second != diagram2[j].second) return (diagram2[i].second < diagram2[j].second); else return (diagram2[i].first > diagram2[j].first); } );
+
+ // Find the inverses of the orders.
+ std::vector<int> order1(num_pts_dgm); std::vector<int> order2(num_pts_dgm);
+ for(int i = 0; i < num_pts_dgm; i++) for (int j = 0; j < num_pts_dgm; j++) if(orderp1[j] == i){ order1[i] = j; break; }
+ for(int i = 0; i < num_pts_dgm; i++) for (int j = 0; j < num_pts_dgm; j++) if(orderp2[j] == i){ order2[i] = j; break; }
+
+ // Record all inversions of points in the orders as theta varies along the positive half-disk.
+ std::vector<std::vector<std::pair<int,double> > > anglePerm1(num_pts_dgm);
+ std::vector<std::vector<std::pair<int,double> > > anglePerm2(num_pts_dgm);
+
+ int m1 = angles1.size();
+ for (int i = 0; i < m1; i++){
+ double theta = angles1[i].first; int p = angles1[i].second.first; int q = angles1[i].second.second;
+ anglePerm1[order1[p]].emplace_back(p,theta);
+ anglePerm1[order1[q]].emplace_back(q,theta);
+ int a = order1[p]; int b = order1[q]; order1[p] = b; order1[q] = a;
+ }
+
+ int m2 = angles2.size();
+ for (int i = 0; i < m2; i++){
+ double theta = angles2[i].first; int p = angles2[i].second.first; int q = angles2[i].second.second;
+ anglePerm2[order2[p]].emplace_back(p,theta);
+ anglePerm2[order2[q]].emplace_back(q,theta);
+ int a = order2[p]; int b = order2[q]; order2[p] = b; order2[q] = a;
+ }
+
+ for (int i = 0; i < num_pts_dgm; i++){
+ anglePerm1[order1[i]].emplace_back(i,pi/2);
+ anglePerm2[order2[i]].emplace_back(i,pi/2);
+ }
+
+ // Compute the SW distance with the list of inversions.
+ for (int i = 0; i < num_pts_dgm; i++){
+ std::vector<std::pair<int,double> > u,v; u = anglePerm1[i]; v = anglePerm2[i];
+ double theta1, theta2; theta1 = -pi/2;
+ unsigned int ku, kv; ku = 0; kv = 0; theta2 = std::min(u[ku].second,v[kv].second);
+ while(theta1 != pi/2){
+ if(diagram1[u[ku].first].first != diagram2[v[kv].first].first || diagram1[u[ku].first].second != diagram2[v[kv].first].second)
+ if(theta1 != theta2)
+ sw += compute_int(theta1, theta2, u[ku].first, v[kv].first, diagram1, diagram2);
+ theta1 = theta2;
+ if ( (theta2 == u[ku].second) && ku < u.size()-1 ) ku++;
+ if ( (theta2 == v[kv].second) && kv < v.size()-1 ) kv++;
+ theta2 = std::min(u[ku].second, v[kv].second);
+ }
+ }
+ }
+
+
+ else{
+ double step = pi/approx;
+
+ // Add projections onto diagonal.
+ int n1, n2; n1 = diagram1.size(); n2 = diagram2.size();
+ for (int i = 0; i < n2; i++)
+ diagram1.emplace_back( (diagram2[i].first + diagram2[i].second)/2, (diagram2[i].first + diagram2[i].second)/2 );
+ for (int i = 0; i < n1; i++)
+ diagram2.emplace_back( (diagram1[i].first + diagram1[i].second)/2, (diagram1[i].first + diagram1[i].second)/2 );
+ int n = diagram1.size();
+
+ // Sort and compare all projections.
+ #ifdef GUDHI_USE_TBB
+ tbb::parallel_for(0, 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) );
+ l2.emplace_back( j, diagram2[j].first*cos(-pi/2+i*step) + diagram2[j].second*sin(-pi/2+i*step) );
+ }
+ std::sort(l1.begin(),l1.end(), [=](const std::pair<int,double> & p1, const std::pair<int,double> & p2){return p1.second < p2.second;});
+ std::sort(l2.begin(),l2.end(), [=](const std::pair<int,double> & p1, const std::pair<int,double> & p2){return p1.second < p2.second;});
+ double f = 0; for (int j = 0; j < n; j++) f += std::abs(l1[j].second - l2[j].second);
+ sw += f*step;
+ });
+ #else
+ for (int i = 0; i < 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) );
+ l2.emplace_back( j, diagram2[j].first*cos(-pi/2+i*step) + diagram2[j].second*sin(-pi/2+i*step) );
+ }
+ std::sort(l1.begin(),l1.end(), [=](const std::pair<int,double> & p1, const std::pair<int,double> & p2){return p1.second < p2.second;});
+ std::sort(l2.begin(),l2.end(), [=](const std::pair<int,double> & p1, const std::pair<int,double> & p2){return p1.second < p2.second;});
+ double f = 0; for (int j = 0; j < n; j++) f += std::abs(l1[j].second - l2[j].second);
+ sw += f*step;
+ }
+ #endif
+ }
+
+ return sw/pi;
+ }
+
+
+ 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 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);
+ }
+
+
+};
+
+} // namespace Sliced_Wasserstein
+} // namespace Gudhi
+
+#endif // SLICED_WASSERSTEIN_H_