QuantumWalkP2

class QuantumWalkP2(N)[source]

Bases: object

Methods providing a very basic framework of a 2 particle CTQW.

Caution

This object constructs some of the framework of a 2 particle quantum walk, but is not feature complete. For a feature complete construction, see pyCTQW.MPI.Graph2P() or pyCTQW.MPI.Line2P() for which this is a base.

Method Summary

QuantumWalkP2.destroy() Destroys the 2 particle quantum walk object, and all associated PETSc matrices/vectors.
QuantumWalkP2.exportPartialTrace(filename, ...) Exports the partial trace of the density matrix to a file.
QuantumWalkP2.exportState(filename, filetype) Exports the final state psi of the quantum walk to a file.
QuantumWalkP2.importInitState(filename, filetype) Imports the initial state of the quantum walk from a file.
QuantumWalkP2.marginal(vec)
QuantumWalkP2.plotEntanglement(filename) Creates a plot of the entanglement over time.
QuantumWalkP2.plotNode(filename, node[, t]) Creates a plot of the marginal probablities on a specified node over time.
QuantumWalkP2.plotNodes(filename[, p, t]) Creates a plot of the node probablities over time.
QuantumWalkP2.propagate(t[, method]) Propagates the quantum walk for time t.
QuantumWalkP2.psiToInit() Copies the state self.psi to self.psi0.
QuantumWalkP2.watch(nodes[, watchtype]) Creates a handle that watches either marginal probability or walker entanglement during propagation.

Method Details

QuantumWalkP2.destroy()[source]

Destroys the 2 particle quantum walk object, and all associated PETSc matrices/vectors.

QuantumWalkP2.exportPartialTrace(filename, filetype, p=1)[source]

Exports the partial trace of the density matrix to a file.

Parameters:
  • filename (str) – path to desired output file.
  • filetype (str) –

    the filetype of the exported adjacency matrix.

    • 'txt' - an \(N\times N\) 2D array in text format.
    • 'bin' - an \(N\times N\) PETSc binary matrix.
  • p (int) – (1|2) - the particle to trace over.

For example, when p=1, \(\rho_2 = \text{Tr}_1(\rho_{12})\), where \(\rho_{12} = |\psi(t)\rangle\langle\psi(t)|\), is exported.

Fortran interface

This function calls the Fortran function partial_trace_mat() for bin filetypes, and the function partial_trace_array() otherwise.

QuantumWalkP2.exportState(filename, filetype)[source]

Exports the final state psi of the quantum walk to a file.

Parameters:
  • filename (str) – path to desired output file.
  • filetype (str) –

    the filetype of the exported adjacency matrix.

    • 'txt' - an \(N\times N\) 2D array in text format.
    • 'bin' - an \(N^2\) element PETSc binary vector.
QuantumWalkP2.importInitState(filename, filetype)[source]

Imports the initial state of the quantum walk from a file.

Parameters:
  • filename (str) – path to the file containing the input state.
  • filetype (str) –

    the filetype of the imported adjacency matrix.

    • 'txt' - an \(N\times N\) 2D array in text format.
    • 'bin' - an \(N^2\) element PETSc binary vector.
Returns:

: this creates a PETSc vector containing the initial statespace, accessed via self.psi0.

Return type:

petsc4py.PETSc.Vec()

Warning

The number of elements in the imported input statespace must equal \(N^2\) (either as a \(N\times N\) text file or a \(N^2\) PETSc binary vector), where \(N\) is the number of nodes the 2P CTQW object is initialized with.

QuantumWalkP2.marginal(vec)[source]
QuantumWalkP2.plotEntanglement(filename)[source]

Creates a plot of the entanglement over time.

Parameters:filename (str) – the absolute/relative path to the desired output file.

Important

watch() must be called prior to propagation in order for the entanglement to be stored during propagation.

Note

  • ensure a file extension is present so that filetype is correctly set (choose one of png, pdf, ps, eps or svg).
  • Timesteps are recorded when propagations occur - to increase the number of timesteps, increase the number of propagations prior to plotting.
QuantumWalkP2.plotNode(filename, node, t=None)[source]

Creates a plot of the marginal probablities on a specified node over time.

Parameters:
  • filename (str) – the absolute/relative path to the desired output file.
  • node (int) – the node to plot.

Important

watch() must be called prior to propagation in order for the probabilities at the specified nodes to be stored.

Note

  • ensure a file extension is present so that filetype is correctly set (choose one of png, pdf, ps, eps or svg).
  • Timesteps are recorded when propagations occur - to increase the number of timesteps, increase the number of propagations.
  • See plotNodes() to plot multiple nodes.
QuantumWalkP2.plotNodes(filename, p=1, t=None)[source]

Creates a plot of the node probablities over time.

Parameters:
  • filename (str) – the absolute/relative path to the desired output file.
  • p (int) – (1|2) - choose whether to plot the marginal probability of particle 1 or 2.

Important

watch() must be called prior to propagation in order for the probabilities at the specified nodes to be stored.

Note

  • if multiple nodes are watched, they will all be plotted.
  • if you wish to plot marginal probabilities for both particle 1 and 2, see plotNode().
  • ensure a file extension is present so that filetype is correctly set (choose one of png, pdf, ps, eps or svg).
  • Timesteps are recorded when propagations occur - to increase the number of timesteps, increase the number of propagations.
QuantumWalkP2.propagate(t, method='chebyshev', **kwargs)[source]

Propagates the quantum walk for time t.

Parameters:
  • t (float) – the timestep over which to propagate the state walk.psi0 via the relationship \(\left|\psi\right\rangle = e^{-iHt}\left|\psi0\right\rangle\).
  • method (str) – the propagation algorithm to use ('chebyshev' (default) or 'krylov').
Keyword Arguments:
 
  • esolver (str) – the eigensolver algorithm to use.

    • 'krylovschur' (default) - Krylov-Schur
    • 'arnoldi' - Arnoldi Method
    • 'lanczos' - Lanczos Method
    • 'power' - Power Iteration/Rayleigh Quotient Iteration
    • 'gd' - Generalized Davidson
    • 'jd' - Jacobi-Davidson,
    • 'lapack' - Uses LAPACK eigenvalue solver routines
    • 'arpack' - only available if SLEPc is compiled with ARPACK linking
  • workType (str) – can be used to set the eigensolver worktype (either 'ncv' or 'mpd'). The default is to let SLEPc decide.

  • workSize (int) – sets the work size if workType is set.

  • tolIn (float) – tolerance of the eigenvalue solver (default 0. (SLEPc decides)).

  • maxIt (int) – maximum number of iterations of the eigenvalue solver (default 0 (SLEPc decides)).

  • verbose (bool) – if True, writes eigensolver information to the console

  • emax_estimate (float) – used to override the calculation of the graphs maximum eigenvalue.

  • emin_estimate (float) – used to override the calculation of the graphs minimum eigenvalue. For a finite graph, emin_estimate=0 is set by default.

    Caution

    • If supplied, the value of emax_estimate (\(\hat{\lambda}_{\max}\)) must satisfy \(\hat{\lambda}_{\max}\geq\lambda_{\max}\), where \(\lambda_{\max}\) is the actual maximum eigenvalue of the graph.
    • Similarly, the value of emin_estimate (\(\hat{\lambda}_{\min}\)) must satisfy \(\hat{\lambda}_{\min}\leq\lambda_{\min}\), where \(\lambda_{\min}\) is the actual minimum eigenvalue of the graph.
    • The greater the difference value of \(|\hat{\lambda}_{\max} -\lambda_{\max}|\) and \(|\hat{\lambda}_{\min} -\lambda_{\min}|\), the longer the convergence time of the chebyshev propagator.

    Note

    • The keyword arguments properties only apply if propagator='chebyshev'
    • For more information regarding these properties,refer to the SLEPc documentation
Returns:

: this creates a PETSc vector containing the propagated state, accessed via the attribute psi, as well as a PETSc vector containing the marginal probabilities, prob.

Return type:

petsc4py.PETSc.Vec()

Note

The EigSolver properties only take effect if method=’chebyshev’.

Fortran interface

This function calls the Fortran function expm(). for the Krylov method, and the function qw_cheby() for the Chebyshev method.

QuantumWalkP2.psiToInit()[source]

Copies the state self.psi to self.psi0.

Example

This can be used to iterate the propagation over discrete timesteps:

for i in range(10):
        walk.propagate(0.01,method='chebyshev')
        walk.psiToInit()
QuantumWalkP2.watch(nodes, watchtype='prob', **kwargs)[source]

Creates a handle that watches either marginal probability or walker entanglement during propagation.

Parameters:
  • nodes (array of ints) – the nodes to watch marginal probability (e.g. [0,1,4]).
  • watchtype (str) –

    (‘prob’ , ‘entanglement’).

    the type of watch handle to produce.

Keyword Arguments:
 

: If ``watchtype=’entanglement’``, EigSolver keywords can also be passed;

for more details of the available EigSolver properties, see propagate().

Returns:

: * if watchtype='prob', creates a handle that can be accessed to retrieve marginal node probabilities for various \(t\) * if watchtype='entanglment', creates a handle that can be accessed to retrieve entanglement values for various \(t\)

Return type:

Examples

To watch the entanglement,
>>> walk.watch(None,watchtype='entanglement')
>>> walk.propagate(5.,method='chebyshev')
>>> timeArray, entArray = walk.entanglementHandle.getEntanglement()

Note

  • The entanglement measure used in Von Neumann entropy, calculated via \(S=-\sum_{i}\lambda_i\log_{2}\lambda_i\), where \(\lambda_i\) are the eigenvalues of the reduced density matrix \(\rho_2 = \text{Tr}_1(|\psi(t)\rangle\langle\psi(t)|)\)
  • Nodes do not need to be specified, as entanglement is a global measurement.
  • As it is a global measurement, there is a large amount of node communication which may increase overall program run time.
To watch the probabilities,
>>> walk.watch([0,1,4])
>>> walk.propagate(2.,method='chebyshev')
>>> timeArray, probXArray, probYArray = walk.handle.getLocalNode(4,p=2)

Warning

Note that walk.handle attributes are not collective; if running on multiple nodes, only local values will be returned.