Title: | Classical Boson Sampling |
---|---|
Description: | Classical Boson Sampling using the algorithm of Clifford and Clifford (2017) <arXiv:1706.01260>. Also provides functions for generating random unitary matrices, evaluation of matrix permanents (both real and complex) and evaluation of complex permanent minors. |
Authors: | Peter Clifford [aut], Raphaël Clifford [cre, aut] |
Maintainer: | Raphaël Clifford <[email protected]> |
License: | GPL-2 |
Version: | 0.1.5 |
Built: | 2024-12-11 07:23:11 UTC |
Source: | CRAN |
Classical Boson Sampling using the algorithm of Clifford and Clifford (2017) <arXiv:1706.01260>. Also provides functions for generating random unitary matrices, evaluation of matrix permanents (both real and complex) and evaluation of complex permanent minors.
Index of help topics:
BosonSampling-package Classical Boson Sampling Permanent-functions Functions for evaluating matrix permanents bosonSampler Function for independently sampling from the Boson Sampling distribution randomUnitary Random unitary
Peter Clifford <[email protected]> and Raphaël Clifford <[email protected]>
Maintainer: Raphaël Clifford <[email protected]>
The function implements the Boson Sampling algorithm defined in Clifford and Clifford (2017) https://arxiv.org/abs/1706.01260
bosonSampler(A, sampleSize, perm = FALSE)
bosonSampler(A, sampleSize, perm = FALSE)
A |
the first |
sampleSize |
the number of independent sample values required for given |
perm |
TRUE if the permanents and pmfs of each sample value are required |
Let the matrix A
be the first n
columns of an (m x m)
random unitary matrix, then
X <- bosonSampler(A, sampleSize = N, perm = TRUE)
provides X$values
, X$perms
and X$pmfs
,
The component X$values
is an (n x N)
matrix with columns that are
independent sample values from the Boson Sampling distribution.
Each sample value is a vector of n
integer-valued output modes in random order. The elements of the vector can be sorted in
increasing order to provide a multiset representation of the sample value.
The outputs X$perms
and X$pmfs
are
vectors of the permanents and probability mass functions (pmfs) associated with the sample values.
The permanent associated with a sample value v = (v_1,...,v_n)
is the permanent of an (n x n)
matrix constructed with rows
v_1,...,v_n
of A
. Note the constructed matrix, M
, may have repeated rows since v_1,...,v_n
are not necessarily distinct.
The pmf is calculated as Mod(pM)^2/prod(factorial(tabulate(c))
where pM
is the permanent of M
.
X = bosonSampler(A, sampleSize = N, perm = TRUE)
provides X$values
, X$perms
and X$pmfs
. See Details.
Clifford, P. and Clifford, R. (2017) The Classical Complexity of Boson Sampling, https://arxiv.org/abs/1706.01260
set.seed(7) n <- 20 # number of photons m <- 200 # number of output modes A <- randomUnitary(m)[,1:n] # sample of output vectors valueList <- bosonSampler(A, sampleSize = 10)$values valueList # sample of output multisets apply(valueList,2, sort) # set.seed(7) n <- 12 # number of photons m <- 30 # number of output modes A <- randomUnitary(m)[,1:n] # sample of output vectors valueList = bosonSampler(A, sampleSize = 1000)$values # Compare frequency of output modes at different # positions in the output vectors matplot(1:m,apply(valueList,1,tabulate), pch =20, t = "p", xlab = "output modes", ylab = "frequency")
set.seed(7) n <- 20 # number of photons m <- 200 # number of output modes A <- randomUnitary(m)[,1:n] # sample of output vectors valueList <- bosonSampler(A, sampleSize = 10)$values valueList # sample of output multisets apply(valueList,2, sort) # set.seed(7) n <- 12 # number of photons m <- 30 # number of output modes A <- randomUnitary(m)[,1:n] # sample of output vectors valueList = bosonSampler(A, sampleSize = 1000)$values # Compare frequency of output modes at different # positions in the output vectors matplot(1:m,apply(valueList,1,tabulate), pch =20, t = "p", xlab = "output modes", ylab = "frequency")
These three functions are used in the classical Boson Sampling problem
cxPerm(A) rePerm(B) cxPermMinors(C)
cxPerm(A) rePerm(B) cxPermMinors(C)
A |
a square complex matrix. |
B |
a square real matrix. |
C |
a rectangular complex matrix where |
Permanents are evaluated using Glynn's formula (equivalently that of Nijenhuis and Wilf (1978))
cxPerm(A)
returns a complex number: the permanent of the complex matrix A
.rePerm(B)
returns a real number: the permanent of the real matrix B
.cxPermMinors(C)
returns a complex vector of length ncol(C)+1
: the permanents of all
ncol(C)
-dimensional square matrices constructed by removing individual rows from C
.
Glynn, D.G. (2010) The permanent of a square matrix. European Journal of Combinatorics, 31(7):1887–1891.
Nijenhuis, A. and Wilf, H. S. (1978). Combinatorial algorithms: for computers and calculators. Academic press.
set.seed(7) n <- 20 A <- randomUnitary(n) cxPerm(A) # B <- Re(A) rePerm(B) # C <- A[,-n] v <- cxPermMinors(C) # # Check Laplace expansion by sub-permanents c(cxPerm(A),sum(v*A[,n]))
set.seed(7) n <- 20 A <- randomUnitary(n) cxPerm(A) # B <- Re(A) rePerm(B) # C <- A[,-n] v <- cxPermMinors(C) # # Check Laplace expansion by sub-permanents c(cxPerm(A),sum(v*A[,n]))
Returns a square complex matrix sampled from the Haar random unitary distribution.
randomUnitary(size)
randomUnitary(size)
size |
dimension of matrix |
A square complex matrix.
m <- 25 # size of matrix (m x m) set.seed(7) U <- randomUnitary(m) # n <- 5 # First n columns A <- U[,1:n]
m <- 25 # size of matrix (m x m) set.seed(7) U <- randomUnitary(m) # n <- 5 # First n columns A <- U[,1:n]