2-Wasserstein distance calculation

Background

The 2-Wasserstein distance \(W\) is a metric to describe the distance between two distributions, representing e.g. two different conditions \(A\) and \(B\).

For continuous distributions, it is given by

\[W := W(F_A, F_B) = \bigg( \int_0^1 \big|F_A^{-1}(u) - F_B^{-1}(u) \big|^2 du \bigg)^\frac{1}{2},\]

where \(F_A\) and \(F_B\) are the corresponding cumulative distribution functions (CDFs) and \(F_A^{-1}\) and \(F_B^{-1}\) the respective quantile functions.

We specifically consider the squared 2-Wasserstein distance \(d := W^2\) which offers the following decomposition into location, size, and shape terms: \[d := d(F_A, F_B) = \int_0^1 \big|F^{-1}(u) - F^{-1}(u) \big|^2 du = \underbrace{\big(\mu_A - \mu_B\big)^2}_{\text{location}} + \underbrace{\big(\sigma_A - \sigma_B\big)^2}_{\text{size}} + \underbrace{2\sigma_A \sigma_B \big(1 - \rho^{A,B}\big)}_{\text{shape}},\]

where \(\mu_A\) and \(\mu_B\) are the respective means, \(\sigma_A\) and \(\sigma_B\) are the respective standard deviations, and \(\rho^{A,B}\) is the Pearson correlation of the points in the quantile-quantile plot of \(F_A\) and \(F_B\).

Usage in two-sample setting

In case the distributions \(F_A\) and \(F_B\) are not explicitly given and information is only availbale in the form of samples from \(F_A\) and \(F_B\), respectively, we use the corresponding empirical CDFs \(\hat{F}_A\) and \(\hat{F}_B\). Then, the 2-Wasserstein distance is computed by

\[d(\hat{F}_A, \hat{F}_B) \approx \frac{1}{K} \sum_{k=1}^K \big(Q_A^{\alpha_k} - Q_B^{\alpha_k} \big) \approx \big(\hat{\mu}_A - \hat{\mu}_B\big)^2 + \big(\hat{\sigma}_A - \hat{\sigma}_B\big)^2 + 2\hat{\sigma}_A \hat{\sigma}_B \big(1 - \hat{\rho}^{A,B}\big).\]

Here, \(Q_A\) and \(Q_B\) denote equidistant quantiles of \(F_A\) and \(F_B\), respectively, at the levels \(\alpha_k := \frac{k-0.5}{K}, k = 1, \dots , K\), where we use \(K:=1000\) in our implementation. Moreover, \(\hat{\mu}_A, \hat{\mu}_B, \hat{\sigma}_A, \hat{\sigma}_B\) and \(\hat{\rho}_{A,B}\) denote the empirical versions of the corresponding quantities.

Three implementations

The waddR package offers three functions to compute the 2-Wasserstein distance in two-sample settings.

We will use samples from normal distributions to illustrate them.

library(waddR)

set.seed(24)
x <- rnorm(100,mean=0,sd=1)
y <- rnorm(100,mean=2,sd=1)

The first function, wasserstein_metric, offers a faster reimplementation in C++ of the wasserstein1d function from the R package transport, which is able to compute general \(p\)-Wasserstein distances. For \(p=2\), we obtain the 2-Wasserstein distance \(W\).

wasserstein_metric(x,y,p=2)
#> [1] 2.044457

The corresponding value of the squared 2-Wasserstein distance \(d\) is then computed as:

wasserstein_metric(x,y,p=2)^2
#> [1] 4.179803

The second function, squared_wass_approx, computes the squared 2-Wasserstein distance by calculating the mean squared difference of the equidistant quantiles (first approximation in the previous formula). This function is currently used to compute the 2-Wasserstein distances in the testing procedures.

squared_wass_approx(x,y)
#> [1] 4.179803

The third function, squared_wass_decomp, approximates the squared 2-Wasserstein distance using the location, size, and shape terms from the above decomposition (second apporximation in the previous formula). It also returns the respective decomposition values.

squared_wass_decomp(x,y)
#> $distance
#> [1] 4.180458
#> 
#> $location
#> [1] 4.114983
#> 
#> $size
#> [1] 0.002307
#> 
#> $shape
#> [1] 0.06316766

In the considered example, the decomposition results suggest that the two distributions differ with respect to location (mean), but not in terms of size and shape, thus confirming the underlying normal model.

See also

Session Info

sessionInfo()
#> R version 4.2.0 RC (2022-04-19 r82224)
#> Platform: x86_64-pc-linux-gnu (64-bit)
#> Running under: Ubuntu 20.04.4 LTS
#> 
#> Matrix products: default
#> BLAS:   /home/biocbuild/bbs-3.15-bioc/R/lib/libRblas.so
#> LAPACK: /home/biocbuild/bbs-3.15-bioc/R/lib/libRlapack.so
#> 
#> locale:
#>  [1] LC_CTYPE=en_US.UTF-8       LC_NUMERIC=C              
#>  [3] LC_TIME=en_GB              LC_COLLATE=C              
#>  [5] LC_MONETARY=en_US.UTF-8    LC_MESSAGES=en_US.UTF-8   
#>  [7] LC_PAPER=en_US.UTF-8       LC_NAME=C                 
#>  [9] LC_ADDRESS=C               LC_TELEPHONE=C            
#> [11] LC_MEASUREMENT=en_US.UTF-8 LC_IDENTIFICATION=C       
#> 
#> attached base packages:
#> [1] stats     graphics  grDevices utils     datasets  methods   base     
#> 
#> other attached packages:
#> [1] waddR_1.10.0
#> 
#> loaded via a namespace (and not attached):
#>  [1] MatrixGenerics_1.8.0        Biobase_2.56.0             
#>  [3] httr_1.4.2                  sass_0.4.1                 
#>  [5] bit64_4.0.5                 jsonlite_1.8.0             
#>  [7] splines_4.2.0               bslib_0.3.1                
#>  [9] assertthat_0.2.1            stats4_4.2.0               
#> [11] BiocFileCache_2.4.0         blob_1.2.3                 
#> [13] GenomeInfoDbData_1.2.8      yaml_2.3.5                 
#> [15] pillar_1.7.0                RSQLite_2.2.12             
#> [17] lattice_0.20-45             glue_1.6.2                 
#> [19] digest_0.6.29               GenomicRanges_1.48.0       
#> [21] XVector_0.36.0              minqa_1.2.4                
#> [23] htmltools_0.5.2             Matrix_1.4-1               
#> [25] pkgconfig_2.0.3             zlibbioc_1.42.0            
#> [27] purrr_0.3.4                 BiocParallel_1.30.0        
#> [29] lme4_1.1-29                 arm_1.12-2                 
#> [31] tibble_3.1.6                generics_0.1.2             
#> [33] IRanges_2.30.0              ellipsis_0.3.2             
#> [35] cachem_1.0.6                SummarizedExperiment_1.26.0
#> [37] BiocGenerics_0.42.0         cli_3.3.0                  
#> [39] magrittr_2.0.3              crayon_1.5.1               
#> [41] memoise_2.0.1               evaluate_0.15              
#> [43] fansi_1.0.3                 nlme_3.1-157               
#> [45] MASS_7.3-57                 tools_4.2.0                
#> [47] lifecycle_1.0.1             matrixStats_0.62.0         
#> [49] stringr_1.4.0               S4Vectors_0.34.0           
#> [51] DelayedArray_0.22.0         eva_0.2.6                  
#> [53] compiler_4.2.0              jquerylib_0.1.4            
#> [55] GenomeInfoDb_1.32.0         rlang_1.0.2                
#> [57] grid_4.2.0                  RCurl_1.98-1.6             
#> [59] nloptr_2.0.0                rappdirs_0.3.3             
#> [61] SingleCellExperiment_1.18.0 bitops_1.0-7               
#> [63] rmarkdown_2.14              boot_1.3-28                
#> [65] abind_1.4-5                 DBI_1.1.2                  
#> [67] curl_4.3.2                  R6_2.5.1                   
#> [69] knitr_1.38                  dplyr_1.0.8                
#> [71] fastmap_1.1.0               bit_4.0.4                  
#> [73] utf8_1.2.2                  filelock_1.0.2             
#> [75] stringi_1.7.6               parallel_4.2.0             
#> [77] Rcpp_1.0.8.3                vctrs_0.4.1                
#> [79] dbplyr_2.1.1                tidyselect_1.1.2           
#> [81] xfun_0.30                   coda_0.19-4