Title: | Handle Continuous-Time Quantum Walks with R |
---|---|
Description: | Functions and tools for creating, visualizing, and investigating properties of continuous-time quantum walks, including efficient calculation of matrices such as the mixing matrix, average mixing matrix, and spectral decomposition of the Hamiltonian. E. Farhi (1997): <arXiv:quant-ph/9706062v2>; C. Godsil (2011) <arXiv:1103.2578v3>. |
Authors: | Vitor Marques [aut, cre, cph] |
Maintainer: | Vitor Marques <[email protected]> |
License: | MIT + file LICENSE |
Version: | 0.1.0.9000 |
Built: | 2024-11-01 03:23:48 UTC |
Source: | https://github.com/vitormarquesr/qwalkr |
Apply a Function to an Operator
act_eigfun(object, ...)
act_eigfun(object, ...)
object |
a representation of the operator. |
... |
further arguments passed to or from other methods. |
The resulting operator from the application of the function.
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) act_eigfun(s, function(x) x^2) #-> act_eigfun.spectral(...)
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) act_eigfun(s, function(x) x^2) #-> act_eigfun.spectral(...)
Apply a function to a Hermitian matrix based on the representation
given by class spectral
.
## S3 method for class 'spectral' act_eigfun(object, FUN, ...)
## S3 method for class 'spectral' act_eigfun(object, FUN, ...)
object |
an instance of class |
FUN |
the function to be applied to the matrix. |
... |
further arguments passed on to |
The matrix resulting from the application of FUN
.
A Hermitian Matrix admits the spectral decomposition
where are its eigenvalues and
the
orthogonal projector onto the
-eigenspace.
If =
FUN
is defined on the eigenvalues of H
, then
act_eigfun
performs the following calculation
H <- matrix(c(0,1,1,1,0,1,1,1,0), nrow=3) decomp <- spectral(H) # Calculates H^2. act_eigfun(decomp, FUN = function(x) x^2) # Calculates sin(H). act_eigfun(decomp, FUN = function(x) sin(x)) # Calculates H^3. act_eigfun(decomp, FUN = function(x, y) x^y, 3)
H <- matrix(c(0,1,1,1,0,1,1,1,0), nrow=3) decomp <- spectral(H) # Calculates H^2. act_eigfun(decomp, FUN = function(x) x^2) # Calculates sin(H). act_eigfun(decomp, FUN = function(x) sin(x)) # Calculates H^3. act_eigfun(decomp, FUN = function(x, y) x^y, 3)
The Average Mixing Matrix of a Quantum Walk
avg_matrix(object, ...)
avg_matrix(object, ...)
object |
a representation of the quantum walk. |
... |
further arguments passed to or from other methods. |
The average mixing matrix.
mixing_matrix()
, gavg_matrix()
, avg_matrix.ctqwalk()
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) avg_matrix(w) #-> avg_matrix.ctqwalk(...)
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) avg_matrix(w) #-> avg_matrix.ctqwalk(...)
The Average Mixing Matrix of a Continuous-Time Quantum Walk
## S3 method for class 'ctqwalk' avg_matrix(object, ...)
## S3 method for class 'ctqwalk' avg_matrix(object, ...)
object |
a representation of the quantum walk. |
... |
further arguments passed to or from other methods. |
Let be the mixing matrix of the quantum walk, then the average mixing matrix is defined as
and encodes the long-term average behavior of the walk. Given the Hamiltonian
, it is possible to prove that
avg_matrix()
returns the average mixing matrix
as a square matrix of the same order as the walk.
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Return the average mixing matrix avg_matrix(walk)
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Return the average mixing matrix avg_matrix(walk)
Returns the adjacency matrix of the cartesian product of two graphs
given the adjacency matrix of each one, and
.
cartesian(G, H = NULL)
cartesian(G, H = NULL)
G |
adjacency matrix of the first graph. |
H |
adjacency matrix of the second graph. If not provided,
it takes the same value as |
Let be the adjacency matrices of
the graphs
such that
and
,
then the adjacency matrix of the cartesian product
is
given by
P3 <- matrix(c(0,1,0,1,0,1,0,1,0), nrow=3) K3 <- matrix(c(0,1,1,1,0,1,1,1,0), nrow=3) # Return the adjacency matrix of P3 X K3 cartesian(P3, K3) # Return the adjacency matrix of P3 X P3 cartesian(P3)
P3 <- matrix(c(0,1,0,1,0,1,0,1,0), nrow=3) K3 <- matrix(c(0,1,1,1,0,1,1,1,0), nrow=3) # Return the adjacency matrix of P3 X K3 cartesian(P3, K3) # Return the adjacency matrix of P3 X P3 cartesian(P3)
ctqwalk()
creates a quantum walk object from a hamiltonian.
ctqwalk(hamiltonian, ...)
ctqwalk(hamiltonian, ...)
hamiltonian |
a Hermitian Matrix representing the Hamiltonian of the system. |
... |
further arguments passed on to |
A list with the walk related objects, i.e the hamiltonian and its spectral
decomposition (See spectral()
for further details)
spectral()
, unitary_matrix.ctqwalk()
,
mixing_matrix.ctqwalk()
, avg_matrix.ctqwalk()
,
gavg_matrix.ctqwalk()
# Creates a walk from the adjacency matrix of the graph P3. ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3))
# Creates a walk from the adjacency matrix of the graph P3. ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3))
The Generalized Average Mixing Matrix of a Quantum Walk
gavg_matrix(object, ...)
gavg_matrix(object, ...)
object |
a representation of the quantum walk. |
... |
further arguments passed to or from other methods. |
The generalized average mixing matrix.
mixing_matrix()
, avg_matrix()
, gavg_matrix.ctqwalk()
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) gavg_matrix(w, rnorm(100)) #-> gavg_matrix.ctqwalk(...)
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) gavg_matrix(w, rnorm(100)) #-> gavg_matrix.ctqwalk(...)
The Generalized Average Mixing Matrix of a Continuous-Time Quantum Walk
## S3 method for class 'ctqwalk' gavg_matrix(object, R, ...)
## S3 method for class 'ctqwalk' gavg_matrix(object, R, ...)
object |
a representation of the quantum walk. |
R |
samples from the random variable |
... |
further arguments passed to or from other methods. |
Let be the mixing matrix of the quantum walk and
a random variable
with associated probability density function
. Then the generalized average mixing
matrix under
is defined as
gavg_matrix()
returns the generalized average mixing matrix
as a square matrix of the same order as the walk.
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Return the average mixing matrix under a Standard Gaussian distribution gavg_matrix(walk, rnorm(1000))
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Return the average mixing matrix under a Standard Gaussian distribution gavg_matrix(walk, rnorm(1000))
Extract an Eigen-Projector from an operator
get_eigproj(object, ...)
get_eigproj(object, ...)
object |
a representation of the operator. |
... |
further arguments passed to or from other methods. |
A representation of the requested eigen-projector.
get_eigspace()
, get_eigschur()
,
get_eigproj.spectral()
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigproj(s, 1) #-> get_eigproj.spectral(...)
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigproj(s, 1) #-> get_eigproj.spectral(...)
Get the orthogonal projector associated with an eigenspace based on the representation
of a Hermitian Matrix given by class spectral
.
## S3 method for class 'spectral' get_eigproj(object, id, ...)
## S3 method for class 'spectral' get_eigproj(object, id, ...)
object |
an instance of class |
id |
index for the desired eigenspace according to the ordered (decreasing) spectra. |
... |
further arguments passed to or from other methods. |
The orthogonal projector of the desired eigenspace.
A Hermitian matrix S
admits the spectral decomposition
such that
is the orthogonal projector onto the
-eigenspace. If
is the matrix associated to the eigenspace, then
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the projector associated to the eigenvalue -1. get_eigproj(decomp, id=2) # Returns the projector associated to the eigenvalue 2. get_eigproj(decomp, id=1)
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the projector associated to the eigenvalue -1. get_eigproj(decomp, id=2) # Returns the projector associated to the eigenvalue 2. get_eigproj(decomp, id=1)
Extract a Schur Cross-Product from an Operator
get_eigschur(object, ...)
get_eigschur(object, ...)
object |
a representation of the operator. |
... |
further arguments passed to or from other methods. |
A representation of the requested Schur cross-product.
get_eigspace()
, get_eigproj()
,
get_eigschur.spectral()
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigschur(s, 1, 2) #-> get_eigschur.spectral(...)
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigschur(s, 1, 2) #-> get_eigschur.spectral(...)
Get the Schur product between eigen-projectors based on the representation of a
Hermitian Matrix given by class spectral
.
## S3 method for class 'spectral' get_eigschur(object, id1, id2 = NULL, ...)
## S3 method for class 'spectral' get_eigschur(object, id1, id2 = NULL, ...)
object |
an instance of class |
id1 |
index for the first eigenspace according to the ordered (decreasing) spectra. |
id2 |
index for the second eigenspace according to the ordered (decreasing) spectra. If not provided,
it takes the same value as |
... |
further arguments passed to or from other methods. |
The Schur product of the corresponding eigenprojectors, .
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the Schur product between the 2-projector and -1-projector. get_eigschur(decomp, id1=2, id2=1) # Returns the Schur square of the 2-projector. get_eigschur(decomp, id1=1, id2=1) # Also returns the Schur square of the 2-projector get_eigschur(decomp, id1=1)
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the Schur product between the 2-projector and -1-projector. get_eigschur(decomp, id1=2, id2=1) # Returns the Schur square of the 2-projector. get_eigschur(decomp, id1=1, id2=1) # Also returns the Schur square of the 2-projector get_eigschur(decomp, id1=1)
Extract an Eigenspace from an Operator
get_eigspace(object, ...)
get_eigspace(object, ...)
object |
a representation of the operator. |
... |
further arguments passed to or from other methods. |
A representation of the requested eigenspace.
get_eigproj()
, get_eigschur()
,
get_eigspace.spectral()
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigspace(s, 1) #-> get_eigspace.spectral(...)
s <- spectral(rbind(c(0.5, 0.3), c(0.3,0.7))) get_eigspace(s, 1) #-> get_eigspace.spectral(...)
Get the eigenbasis associated with an eigenvalue based on the representation
of a Hermitian Matrix given by class spectral
.
## S3 method for class 'spectral' get_eigspace(object, id, ...)
## S3 method for class 'spectral' get_eigspace(object, id, ...)
object |
an instance of class |
id |
index for the desired eigenspace according to the ordered (decreasing) spectra. |
... |
further arguments passed to or from other methods. |
A matrix whose columns form the orthonormal eigenbasis.
If s <- spectral(A)
and V <- s$eigvectors
, then the extracted eigenspace
is some submatrix
V[, _]
.
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the two orthonormal eigenvectors corresponding to the eigenvalue -1. get_eigspace(decomp, id=2) # Returns the eigenvector corresponding to the eigenvalue 2. get_eigspace(decomp, id=1)
# Spectra is {2, -1} with multiplicities one and two respectively. decomp <- spectral(matrix(c(0,1,1,1,0,1,1,1,0), nrow=3)) # Returns the two orthonormal eigenvectors corresponding to the eigenvalue -1. get_eigspace(decomp, id=2) # Returns the eigenvector corresponding to the eigenvalue 2. get_eigspace(decomp, id=1)
Returns the all-ones matrix of order n
.
J(n)
J(n)
n |
the order of the matrix. |
A square matrix of order in which every entry is
equal to 1. The all-ones matrix is given by
.
# Return the all-ones matrix of order 5. J(5)
# Return the all-ones matrix of order 5. J(5)
The Mixing Matrix of a Quantum Walk
mixing_matrix(object, ...)
mixing_matrix(object, ...)
object |
a representation of the quantum walk. |
... |
further arguments passed to or from other methods. |
The mixing matrix of the quantum walk.
unitary_matrix()
, avg_matrix()
, gavg_matrix()
,
mixing_matrix.ctqwalk()
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) mixing_matrix(w, t = 2*pi) #-> mixing_matrix.ctqwalk(...)
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) mixing_matrix(w, t = 2*pi) #-> mixing_matrix.ctqwalk(...)
The Mixing Matrix of a Continuous-Time Quantum Walk
## S3 method for class 'ctqwalk' mixing_matrix(object, t, ...)
## S3 method for class 'ctqwalk' mixing_matrix(object, t, ...)
object |
an instance of class |
t |
it will be returned the mixing matrix at time |
... |
further arguments passed to or from other methods. |
Let be the time evolution operator of the quantum walk at
time
, then the mixing matrix is given by
is a doubly stochastic real symmetric matrix, which encodes the
probability density of the quantum system at time
.
More precisely, the entry gives us the probability
of measuring the standard basis state
at time
, given that
the quantum walk started at
.
mixing_matrix()
returns the mixing matrix of the CTQW
evaluated at time t
.
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Returns the mixing matrix at time t = 2*pi, M(2pi) mixing_matrix(walk, t = 2*pi)
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Returns the mixing matrix at time t = 2*pi, M(2pi) mixing_matrix(walk, t = 2*pi)
Print the ctqwalk output
## S3 method for class 'ctqwalk' print(x, ...)
## S3 method for class 'ctqwalk' print(x, ...)
x |
an object of the class ctqwalk. |
... |
further arguments passed to or from other methods. |
Called mainly for its side effects. However, also
returns x
invisibly.
spectral()
is a wrapper around base::eigen()
designed for Hermitian matrices,
which can handle repeated eigenvalues.
spectral(S, multiplicity = TRUE, tol = .Machine$double.eps^0.5, ...)
spectral(S, multiplicity = TRUE, tol = .Machine$double.eps^0.5, ...)
S |
a Hermitian matrix. Obs: The matrix is always assumed to be Hermitian, and only its lower triangle (diagonal included) is used. |
multiplicity |
if |
tol |
two eigenvalues |
... |
further arguments passed on to |
The spectral decomposition of S
is returned as a list with components
eigvals |
vector containing the unique eigenvalues of |
multiplicity |
multiplicities of the eigenvalues in |
eigvectors |
a |
The Spectral Theorem ensures the eigenvalues of S
are real and that the vector space
admits an orthonormal basis consisting of eigenvectors of S
. Thus, if s <- spectral(S)
,
and V <- s$eigvectors; lam <- s$eigvals
, then
where \Lambda =\
diag(rep(lam, times=s$multiplicity))
base::eigen()
, get_eigspace.spectral()
,
get_eigproj.spectral()
, get_eigschur.spectral()
,
act_eigfun.spectral()
spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Use "tol" to set the tolerance for numerical equality spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3), tol=10e-5) # Use "multiplicity=FALSE" to force each eigenvalue to be considered unique spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3), multiplicity = FALSE)
spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Use "tol" to set the tolerance for numerical equality spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3), tol=10e-5) # Use "multiplicity=FALSE" to force each eigenvalue to be considered unique spectral(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3), multiplicity = FALSE)
Computes the trace of a matrix A
.
tr(A)
tr(A)
A |
a square matrix. |
If has order
,
then
.
A <- rbind(1:5, 2:6, 3:7) # Calculate the trace of A tr(A)
A <- rbind(1:5, 2:6, 3:7) # Calculate the trace of A tr(A)
Computes the trace inner product of two matrices A
and B
.
trdot(A, B)
trdot(A, B)
A , B
|
square matrices. |
The trace inner product on is
defined as
A <- rbind(1:5, 2:6, 3:7) B <- rbind(7:11, 8:12, 9:13) # Compute the trace inner product of A and B trdot(A, B)
A <- rbind(1:5, 2:6, 3:7) B <- rbind(7:11, 8:12, 9:13) # Compute the trace inner product of A and B trdot(A, B)
The Unitary Time Evolution Operator of a Quantum Walk
unitary_matrix(object, ...)
unitary_matrix(object, ...)
object |
a representation of the quantum walk. |
... |
further arguments passed to or from other methods. |
The unitary time evolution operator.
mixing_matrix()
, unitary_matrix.ctqwalk()
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) unitary_matrix(w, t = 2*pi) #-> unitary_matrix.ctqwalk(...)
w <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) unitary_matrix(w, t = 2*pi) #-> unitary_matrix.ctqwalk(...)
The Unitary Time Evolution Operator of a Continuous-Time Quantum Walk
## S3 method for class 'ctqwalk' unitary_matrix(object, t, ...)
## S3 method for class 'ctqwalk' unitary_matrix(object, t, ...)
object |
an instance of class |
t |
it will be returned the evolution operator at time |
... |
further arguments passed to or from other methods. |
If is the quantum state of the system at time
, and
the Hamiltonian operator, then the evolution is governed by
the Schrodinger equation
and if is time-independent its solution is given by
The evolution operator is the result of the complex matrix exponential and it can be calculated as
in which .
unitary_matrix()
returns the unitary time evolution operator of the
CTQW evaluated at time t
.
ctqwalk()
, unitary_matrix()
,
act_eigfun()
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Returns the operator at time t = 2*pi, U(2pi) unitary_matrix(walk, t = 2*pi)
walk <- ctqwalk(matrix(c(0,1,0,1,0,1,0,1,0), nrow=3)) # Returns the operator at time t = 2*pi, U(2pi) unitary_matrix(walk, t = 2*pi)