pyCTQW Manual

Contents Features

Downloads

Latest release

The latest stable version of pyCTQW.MPI is available to download from GitHub; older releases are also available for download and can be found here.

Alternatively, if pip is installed, pyCTQW can be installed using pip:

$ pip install pyCTQW

Please see the installation page for more information on installing pyCTQW.

Development version

If you would like to try out the development version, ths can be downloaded as an archive (pyctqw-master.tar.gz) or alternatively, checked out via Git:

$ git clone https://github.com/josh146/pyctqw.git

or Subversion:

$ svn co https://github.com/josh146/pyctqw

or installed using pip:

$ pip install --allow-unverified pyCTQW pyCTQW==dev

Caution

The development version may be unstable, and can possibly break various functionality from time to time.

Installation

Note

At the moment, only Unix based systems are supported.

Dependencies

In addition to an MPI implementation (e.g. MPICH or Open MPI), the pyCTQW.MPI module depends on the following components:

Installation using pip

After ensuring NumPy and petsc4py are installed (and all PETSc, SLEPc and MPI environment variables are properly set), pyctqw can be installed using pip:

$ pip install pyCTQW

Note

The current development version pyCTQW-dev can also be installed using pip:

$ pip install --allow-unverified pyCTQW pyCTQW==dev

Installing pyCTQW.MPI from source code

Alternatively, the source code can be downloaded and compiled manually:

  1. Download the latest version of pyctqw, extract the pyCTQW archive, and cd into the extracted directory:

    $ wget https://github.com/josh146/pyctqw/archive/1.0.0.tar.gz -O pyctqw-1.0.0.tar.gz
    $ tar xzvf pyctqw-1.0.0.tar.gz
    $ cd pyctqw-1.0.0
    
  2. Ensure that your PETSc and SLEPc environment variables are correctly set; for example,

    $ export PETSC_DIR=/path/to/petsc
    $ export PETSC_ARCH=linux-gnu
    $ export SLEPC_DIR=/path/to/slepc
    

    If you are unsure what your PETSc or SLEPc variables should be, please refer to their documentation.

    Important

    If you plan to install pyCTQW.MPI using root to a system directory, the PETSc and SLEPc environment variables must be available to the root user.

  3. Compile the Python module pyCTQW.MPI by running

    $ python setup.py build
    
  4. System-wide install:

    $ sudo -E python setup.py install
    

    where the command -E ensures that the environment variables set in step 3 are passed to the root.

    Note

    If you do not have root access, or the above command does not appear to work, you can install the package locally by running

    $ python setup.py install --user
    

    Now, have a go running some of the Examples!

Optional: build documentation

If Sphinx is installed, the documentation can be compiled by running

$ make docs

Known Issues

  • Non-mpi fallback modes not present yet

Quickstart Guide

Creating Python Scripts

Once pyCTQW.MPI is installed, calculations involving CTQWs can be called in several ways:

The syntax in both cases are identical, however this guide will be deal with executable Python scripts.


Initialization

The first step is to initialize the PETSc environment, and import the pyCTQW.MPI module:

import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import pyCTQW.MPI as qw

Command Line Options

Next, PETSc can be used to create command line options for the script; for example, the following code creates two command line options, t and N, with default values of 100 and 20 respectively:

OptDB = PETSc.Options()
N = OptDB.getInt('N', 100)
t = OptDB.getReal('t', 20)

Furthermore, most PETSc and SLEPc subroutines accept command line options which modify their settings; for instance, when using the SLEPc EPS eigensolver, the EPS type can be changed dynamically when run through the following options:

$ mpirun -np 2 <program> -eps_type='lapack'

See also

For more details, please refer to the PETSc and SLEPc documentation.


Rank and I/O

When running on multiple nodes, we sometimes want only one node to perform a calculation or operation – for instance, for I/O operations where all nodes already have the same information. Using PETSc, the rank (MPI process number) can be determined for each process, and conditional statements used to control which node performs the I/O operation:

rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
        print '1P Line\n'

Caution

  • Be careful though: all of the methods and functions available in pyCTQW.MPI are designed to work globally on all processes, and should not be created or called on a subset of all available nodes. Doing so may cause your program to hang.
  • Most objects in pyCTQW.MPI contain I/O methods, for instance pyCTQW.MPI.Graph.exportState(); these are global over all nodes (as mentioned above) and should be used over custom methods when possible.
  • I/O operations involving PETSc objects (for instance PETSc.Vec.view()) are also designed to be global/collective over all nodes; see petsc4py for relevant documentation.

Creating CTQW objects

All available objects and methods are detailed over at the API documentation pages pyCTQW.MPI.

See also

You can also refer to the examples to see how CTQW objects can be manipulated


Running your script

Once your script is complete, save it with a .py extension, and make it executable by running

$ chmod +x <script>.py

where <script> is the file path of your script. Then, to run your program, simply run in a terminal

$ mpirun -np X <script>.py [options]

where X is the number of MPI processes to run.


Code profiling

PETSc also allows for easy code profiling by supplying the command line option -log_summary when executing your script, and this is built in to pyCTQW.MPI; for instance, the methods available in pyCTQW.MPI.Line2P() automatically create log stages for creating the Hamiltonian, initial state, finding the eigenvalues, propagation etc.

If you wish to create custom log stages, this can also be done:

stage1 = _PETSc.Log.Stage('First Stage')
stage1.push()
# place stage 1 functions/operations here
stage1.pop()

See also

For more details, please refer to the petsc4py documentation.

pyCTQW.MPI Package

Contains distributed memory methods for calculating:
  • the continuous-time quantum walk on arbitrary graphs and lines for 1, 2 and 3 particles
  • Non-interacting and interacting walks
  • Entanglement calculations
  • ability to place diagonal defects/barriers on graph vertices
  • contains MPI graph isomorphism methods, for both 2 and 3 interacting particles

These classes also contain plotting methods, to visualise the quantum walking dynamics

Arbitrary graph QW

Graph(N[, filename, filetype, d, amp, ...]) Performs and analyses 1 particle continuous-time quantum walks on graphs
Graph2P(N[, filename, filetype, d, amp, ...]) Performs and analyses 2 particle continuous-time quantum walks on graphs
Graph3P(N[, filename, filetype, d, amp, ...]) Performs and analyses 3 particle continuous-time quantum walks on graphs

See also

Examples using these objects and methods:
1P_3cayley.py, 2P_3cayley.py, 3P_3cayley.py

Infinite line QW

Line(N) Performs and analyses 1 particle continuous-time quantum walks on an infinite line
Line2P(N[, d, amp, interaction]) Performs and analyses 2 particle continuous-time quantum walks on an infinite line
Line3P(N[, d, amp, interaction]) Performs and analyses 3 particle continuous-time quantum walks on an infinite line

See also

Examples using these objects and methods:
1P_line.py, 2P_line.py, 3P_line.py

Graph Isomorphism

GraphISO(**kwargs) A graph isomorphism solver, containing functions for creating graph certificates and checking isomorphism of adjacency matrices.

See also

Examples using this object:
GraphISO_sr.py, GraphISO_3cayley.py

Submodules

Caution

The methods and classes available in these submodules have been designed to work with the classes listed above; that is, they are being called implicitly by the above classes, invisible to the user.

It is not recommended that they be used by themselves, however it is possible to do so.

pyCTQW.MPI.io Input and output functions.
pyCTQW.MPI.plot
pyCTQW.MPI.ctqw Additional ctqw related classes.

Fortran Library

Function documentation

libctqwMPI API

(source: ctqwMPI.F90)

I/O subroutines

subroutine exportvec(vec, filename, filetype)

export a PETSc vector to a file

Parameters:
  • vec [real] :: the PETSc vector to export
  • filename [character,in] :: filename only, no path, of the output file
  • filetype [character,in] :: filetype to export (‘txt’ or ‘bin’)
subroutine importvec(vec, filename, filetype)

import a PETSc vector from a file

Parameters:
  • vec [real] :: the PETSc vector to import the file into
  • filename [character,in] :: filename only, no path, of the input file
  • filetype [character,in] :: filetype to export (‘txt’ or ‘bin’)
subroutine exportmat(mat, filename, filetype)

export a PETSc matrix to a file

Parameters:
  • mat [integer] :: the input PETSc matrix
  • filename [character,in] :: the filename of the exported file
  • filetype [character,in] :: filetype to export (‘txt’ or ‘bin’)
subroutine importmat(mat, filename, filetype)

import a PETSc matrix from a file

Parameters:
  • mat [integer] :: the PETSc matrix to import the file into
  • filename [character,in] :: filename only, no path, of the input file
  • filetype [character,in] :: filetype to export (‘txt’ or ‘bin’)

Matrix Array subroutines

subroutine identity(mat, n)

creates a \(n\times n\) identity matrix

Parameters:
  • mat [integer]
  • n [integer] :: the \(n\times n\) identity matrix output
Called from:

adjtoh()

subroutine kron(m1, r1, c1, m2, r2, c2, kronprod)

compute the Kronecker product of two arrays

Parameters:
  • m1 [integer]
  • r1 [real] :: number of rows in matrix 1
  • c1 [real] :: matrix 1
  • m2 [integer]
  • r2 [real] :: number of rows in matrix 2
  • c2 [real] :: kronecker product output
  • kronprod [integer]
Called from:

adjtoh()

Hamiltonian subroutines

subroutine importadjtoh(mat, filename, p, d, amp, interaction, nd)

Import an adjacency matrix from a file, and create a PETSc Hamiltonian matrix.

Parameters:
  • mat [integer] :: output Hamiltonian matrix
  • filename [character,in] :: dense adjacency matrix file in text format
  • p [character,in] :: number of particles in the system (‘1’, ‘2’, or ‘3’)
  • d [real]
  • amp [real]
  • interaction [integer] :: interaction amplitude
  • nd [integer] :: amplitude of defects
Call to:

adjtoh()

subroutine adjtoh(mat, adjarray, p, d, amp, interaction, nd, n)

convert an adjacency array to a PETSc Hamiltonian matrix

Parameters:
  • mat [integer] :: output Hamiltonian matrix
  • adjarray [real]
  • p [character,in] :: number of particles in the system (‘1’, ‘2’, or ‘3’)
  • d [real]
  • amp [real]
  • interaction [integer] :: interaction amplitude
  • nd [integer] :: amplitude of defects
  • n [integer]
Called from:

importadjtoh(), graphiscert(), hamiltonian_p3_line()

Call to:

identity(), kron()

subroutine hamiltonian_p1_line(a, d, amp, nd, n)
Parameters:
  • a [real]
  • d [real]
  • amp [real]
  • nd [integer]
  • n [integer]
subroutine hamiltonian_p2_line(h2, d, amp, interaction, nd, n)
Parameters:
  • h2 [real]
  • d [real]
  • amp [real]
  • interaction [integer]
  • nd [integer]
  • n [integer]
subroutine hamiltonian_p3_line(h3, d, amp, interaction, nd, n)
Parameters:
  • h3 [real]
  • d [real]
  • amp [real]
  • interaction [integer]
  • nd [integer]
  • n [integer]
Call to:

adjtoh()

Statespace subroutines

function coord(x, y, n)

convert from a 2D to 1D statespace for 2 particles, using \(coord = n(x + n/2 - 1) + y + n/2 - 1\)

Parameters:
  • x [integer,in] :: vertex location of particle 1
  • y [integer,in] :: vertex location of particle 2
  • n [integer,in] :: number of nodes in the system
Return:

coord [integer] :: output 2P statepace coordinate

Called from:

p2_init()

function coord3p(x, y, z, n)

convert from a 2D to 1D statespace for 2 particles, using \(coord2P = xn^2 + ny + z\)

Parameters:
  • x [integer,in] :: vertex location of particle 1
  • y [integer,in] :: vertex location of particle 2
  • z [integer,in] :: vertex location of particle 3
  • n [integer,in] :: number of nodes in the system
Return:

coord3p [integer]

Called from:

p3_init()

subroutine marginal1(psi, prob, n)

Calculates the marginal probability of 1 particle state psi; i.e. \(|\psi|^2\)

Parameters:
  • psi [real] :: input statespace PETSc vector
  • prob [real] :: output PETSc vector containing the probabilities
  • n [integer] :: length of vector psi
subroutine marginal2(psi, psim, p, n)

Calculates the marginal probability of particle number p

Parameters:
  • psi [real] :: input PETSc statespace vector of length \(n^p\)
  • psim [real]
  • p [character,in] :: the particle to calculate the marginal probability (‘1’, ‘2’, ‘3’)
  • n [integer] :: number of vertices in the graph
subroutine marginal3(psi, psim, p, n)
Parameters:
  • psi [real]
  • psim [real]
  • p [character,in]
  • n [integer]
subroutine p1_init(psi0, init_state, num, n)
Parameters:
  • psi0 [real]
  • init_state [integer]
  • num [integer]
  • n [integer]
subroutine p2_init(psi0, init_state, num, n)
Parameters:
  • psi0 [real]
  • init_state [integer]
  • num [integer]
  • n [integer]
Called from:

graphiscert()

Call to:

coord()

subroutine p3_init(psi0, init_state, num, n)
Parameters:
  • psi0 [real]
  • init_state [integer]
  • num [integer]
  • n [integer]
Called from:

graphiscert()

Call to:

coord3p()

MatrixExp and Eigenvalues

subroutine qw_krylov(a, t, v, y)
Parameters:
  • a [real]
  • t [real]
  • v [real]
  • y [real]
Called from:

graphiscert()

subroutine min_max_eigs(a, rank_bn, eval, eval_error, which, eig_solver, worktype, worktypeint, tolin, max_it, verbose, error)
Parameters:
  • a [real]
  • rank_bn [real]
  • eval [real]
  • eval_error [real]
  • which [character,in]
  • eig_solver [character,in]
  • worktype [character,in]
  • worktypeint [real]
  • tolin [real]
  • max_it [integer]
  • verbose [real]
  • error [real]
Called from:

graphiscert()

subroutine qw_cheby(psi0, psi, dt, h, emin, emax, rank_bn, n)
Parameters:
  • psi0 [real]
  • psi [real]
  • dt [real]
  • h [real]
  • emin [real]
  • emax [real]
  • rank_bn [real]
  • n [integer]
Called from:

graphiscert()

Entanglement subroutines

subroutine partial_trace_array(psi, rhox, n)
Parameters:
  • psi [real]
  • rhox [real]
  • n [integer]
subroutine partial_trace_mat(psi, rhox, n)
Parameters:
  • psi [real]
  • rhox [real]
  • n [integer]
Called from:

entanglement()

subroutine entanglement(psi, n, vne, eig_solver, worktype, worktypeint, tolin, max_it, verbose, error)
Parameters:
  • psi [real]
  • n [integer]
  • vne [real]
  • eig_solver [character,in]
  • worktype [character,in]
  • worktypeint [real]
  • tolin [real]
  • max_it [integer]
  • verbose [real]
  • error [real]
Call to:

partial_trace_mat()

GraphIso subroutines

function number_of_edges(adjarray, n)
Parameters:
  • adjarray [real]
  • n [integer]
Return:

number_of_edges [integer]

Called from:

getalledgestates(), getalledgestates3p()

function getedgestate(edgenum, adjarray, n)
Parameters:
  • edgenum [real]
  • adjarray [real]
  • n [integer]
Return:

getedgestate [real]

Called from:

getedgestate(), getalledgestates(), getalledgestates3p()

Call to:

getedgestate()

subroutine getalledgestates(init_states, localstatenum, adjarray, n)
Parameters:
  • init_states [integer]
  • localstatenum [integer]
  • adjarray [real]
  • n [integer]
Called from:

graphiscert()

Call to:

number_of_edges(), getedgestate()

subroutine getalledgestates3p(init_states, localstatenum, adjarray, n)
Parameters:
  • init_states [integer]
  • localstatenum [integer]
  • adjarray [real]
  • n [integer]
Called from:

graphiscert()

Call to:

number_of_edges(), getedgestate()

subroutine graphiscert(cert, certlength, adjarray, p, tol, expm_method, eig_solver, emax_estimate, worktype, worktypeint, tolin, max_it, verbose, n)
Parameters:
  • cert [real]
  • certlength [real]
  • adjarray [real]
  • p [real]
  • tol [real]
  • expm_method [character,in]
  • eig_solver [character,in]
  • emax_estimate [real]
  • worktype [character,in]
  • worktypeint [real]
  • tolin [real]
  • max_it [integer]
  • verbose [real]
  • n [integer]
Call to:

adjtoh(), getalledgestates3p(), getalledgestates(), min_max_eigs(), p3_init(), p2_init(), qw_cheby(), qw_krylov(), d_refsor()

subroutine d_refsor(xdont)

:p real xdont(:) [inout]:

Called from:graphiscert()

Note

d_refsor(), a highly optimised Fortran sorting implementation written by Michel Olagnon and part of the ORDERPACK 2.0 suite of ranking and sorting algorithms for Fortran 90.

Caution

It is slightly more difficult to compile and link programs against the Fortran library directly, due to the myriad of pre-processing and the number of required headers/environment variables. As such, it is recommended users interface with this library via the much easier to use Python module pyCTQW.MPI.

However, if you have experience using PETSc and SLEPc with Fortran, and wish to proceed, a makefile template is provided that can be modified to link your program against libctqwMPI.

Compiling libctqwMPI

In addition to an MPI implementation (e.g. MPICH or Open MPI), the Fortran library libctqwMPI depends on the following components:

Once these dependencies are installed, simply open a terminal in the root directory of pyCTQW-X.Y and run

$ make fortran [options]

where available options include

Option Values Description
shared_lib 0 (default), 1

whether to build libctqwMPI as a shared library (shared_lib=1, producing libctqwMPI.so) or a static library (shared_lib=0 (default), producing libctqwMPI.a).

If built as a shared library, compiled programs will be smaller, but libctqwMPI.so will need to be added to a directory used by ld (either by setting the environment variable LD_LIBRARY_PATH or by placing libctqwMPI.so in /usr/local/lib etc).

The fortran library (libctqwMPI.so or libctqwMPI.a) can be found in the pyCTQW-X.Y/lib directory, with required module files found in the pyCTQW-X.Y/include directory.


Using libctqwMPI

To call functions and subroutines from libctqwMPI, your program should have the following structure:

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
program main
    ! interfaces for libctqwMPI
    use ctqwMPI

! PETSc headers required by the preprocessor.
! Note the lack of indentation.
#include <finclude/petsc.h>

    PetscErrorCode :: ierr
    PetscMPIInt    :: rank
    !--------------------
    ! your variable
    ! declarations
    ! -------------------

    ! initialize SLEPc and PETSc
    call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
    call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)

    !--------------------
    ! your program code
    ! -------------------

    ! finalise PETSc
    call PetscFinalize(ierr)

end program main

CTQW subroutines can now be called directly; for example,

call qw_cheby(psi0,psi,t,H,0.*PETSc_i,Emax,rank,n)

See also


Command Line Options

Similarly to the Python scripts, PETSc allows command line arguments to be easily added to your Fortran programs, through the use of the following subroutines:

call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-n",n,flag,ierr)
CHKERRQ(ierr)
if (flag .eqv. .false.) n = 100

In this case, the command line option -n is added, and the command stored in variable n (of type PETSC_INT). Furthermore, if the flag is not provided, the program assigns n a default value of 100.

Furthermore, most PETSc and SLEPc subroutines accept command line options which modify their settings; for instance, when using the SLEPc EPS eigensolver, the EPS type can be changed dynamically when run through the following options:

$ mpirun -np 2 <program> -eps_type='lapack'

See also

For more details, please refer to the PETSc and SLEPc documentation.


Code Profiling

PETSc also allows for easy code profiling using PetscLogStage variables; however, unlike the built in log stages found in pyCTQW.MPI, these must be manually created for fortran programs.

For example, to create enable profiling of they Chebyshev method used for CTQW propagation, we would declare a log stage variable stage,

PetscLogStage  :: stage
PetscLogDouble :: tc0, tc1
PetscErrorCode :: ierr

and then wrap the Chebyshev subroutine like so:

! measure the initial wall time
call PetscTime(tc0,ierr)

! register the log stage, named 'Chebyshev'
call PetscLogStageRegister('Chebyshev',stage,ierr)
call PetscLogStagePush(stage,ierr)

! CTQW propagation
call qw_cheby(psi0,psi,t,H,0.*PETSc_i,Emax,rank,n)

! view the final statespace
call VecView(psi,PETSC_VIEWER_STDOUT_WORLD,ierr)

! wait here until all MPI nodes arrive
call PetscBarrier(psi,ierr)
! end the log stage
call PetscLogStagePop(ierr)

! measure the end wall time
call PetscTime(tc1,ierr)

Then, when calling the program with the -log_summary command line option,

$ mpirun -np 2 <program> -log_summary

the produced code profile will include a summary of performance within the ‘Chebyshev’ log stage.

See also

For more details, please refer to the PETSc and SLEPc documentation.


Linking libctqwMPI

To compile your program and link it against libctqwMPI, you will need to include the ctqw_common file (present in the root directory of the pyCTQW-X.Y folder) in your makefile - this helps set up the required library/include variables.

For example, consider the makefile below (used to compile exampleMPI.F90):

Makefile

[Download makefile]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
#!/usr/bin/make -f

# If they are not automatically found, make sure
# your PETSc and SLEPc variables are correctly
# set by uncommenting out the below.
# are correctly
#PETSC_DIR = 
#PETSC_ARCH =
#SLEPC_DIR = 

CTQW_DIR = ../..
include $(CTQW_DIR)/ctqw_common

program = exampleMPI

#Makefile	
all: $(program)

$(program): $(program).o
	$(FLINKER) $(program).o -o $(program) $(CTQW_LIBS)

$(program).o: $(program).F90
	$(FLINKER) $(program).F90 $(CTQW_INCLUDE_DIRS) -c -o $(program).o 

clean::
	rm -f *.o *.mod *.s

.DEFAULT_GOAL := all

A brief summary of the variables:

CTQW_DIR the path to the pyCTQW-X.Y directory (used to correctly set up CTQW_LIBS and CTQW_INCLUDE_DIRS).
CTQW_LIBS

the libraries required to link against libctqwMPI; these include libctqwMPI.so (by default, this is compiled to CTQW_DIR/lib) as well as the PETSc and SLEPc libraries.

These are required when linking your object files.

CTQW_INCLUDE_DIRS

the location of the header/module files required - by default, these are compiled to CTQW_DIR/include. Note that this also contains required PETSc/SLEPc include directories.

These are required when compiling your code to an object file.


Running your program

To run your program, simply run in a terminal

$ mpirun -np X <program> [options]

where X is the number of MPI nodes you would like to run it on.

Important

If you compiled your program using the shared library libctqwMPI.so, make sure that the shared library is available at runtime, either by:

  • setting the environment variable LD_LIBRARY_PATH to include the directory containing libctqwMPI.so,
  • or by placing libctqwMPI.so in a directory used by ld (e.g. /usr/local/lib etc).

Acknowledgements

The graph isomorphism subroutine graphiscert() uses the external subroutine d_refsor(), a highly optimised Fortran sorting implementation written by Michel Olagnon and part of the ORDERPACK 2.0 suite of ranking and sorting algorithms for Fortran 90.

Examples

CTQW Examples

1P_line.py

Description

This example propagates a 1 particle continuous-time quantum walk on an infinite line

Amongst the features used, it illustrates:
  • recieving command line options using PETSc

  • the use of the chebyshev algorithm
    • setting the EigSolver tolerance
  • adding a diagonal defects to various nodes

  • creating node handles to watch the probability at specified nodes

  • various plotting abilities:
    • probability vs node plots
    • probability vs time plots
Output
_images/1p_line_plot.png _images/1p_line_nodes.png
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 100)
t = OptDB.getReal('t', 20)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '1P Line\n'

# initialise an N (default 100) node graph CTQW
walk = qw.Line(N)

# Create a Hamiltonian with defect and amplitude as below.
d = [3,4]
amp = [2.0,1.5]
walk.createH(d,amp)

# create the initial state (1/sqrt(2)) (|0>+|1>)
init_state = [[0.,1.0/np.sqrt(2.0)], [1.,1.0/np.sqrt(2.0)]]
walk.createInitState(init_state)

# set the eigensolver properties.
walk.EigSolver.setEigSolver(tol=1.e-3)

# create a handle to watch the probability at nodes -5,0,1:
walk.watch([0,1,-5])

# Propagate the CTQW using the Chebyshev method
# for t=100s in timesteps of dt=0.01
# Note that psiToInit() is being used rather than global timesteps.
for i in range(int(t/0.01)):
	walk.propagate(0.01,method='chebyshev')
	walk.psiToInit()

# plot the marginal probabilities
# after propagation over all nodes
walk.plot('out/1p_line_plot.png')

# plot the probability over time for the watched nodes
walk.plotNodes('out/1p_line_nodes.png')

# export final state
walk.exportState("out/1p_final_state.txt", "txt")

# destroy the quantum walk
walk.destroy()

2P_line.py

Description

This example propagates a 2 particle continuous-time quantum walk on an infinite line

Amongst the features used, it illustrates:
  • recieving command line options using PETSc

  • the use of the chebyshev algorithm
    • setting the EigSolver tolerance, as well as the minimum eigenvalue
  • adding a diagonal defects to various nodes

  • same-node interactions between particles

  • creating node handles to watch the probability at specified nodes

  • creating entanglement handles to watch the entanglement

  • various plotting abilities:
    • probability vs node plots
    • probability vs time plots
    • entanglement vs time plot
  • Exporting the final state to a PETSc binary vector file

Output
_images/2p_line_plot.png _images/2p_line_nodes.png _images/2p_line_ent.png
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 100)
t = OptDB.getReal('t', 20)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '2P Line\n'

# initialise an N (default 100) node graph CTQW
walk = qw.Line2P(N)

# Create a Hamiltonian with 2P interaction.
walk.createH(interaction=1.)

# create the initial state (1/sqrt(2)) (|0,0>-|1,1>)
init_state = [[0,0,1.0/np.sqrt(2.0)], [1,1,-1.0/np.sqrt(2.0)]]
walk.createInitState(init_state)

# set the eigensolver properties.
walk.EigSolver.setEigSolver(tol=1.e-2)
# underestimate the minimum eigenvalue
walk.EigSolver.setEigSolver(emin_estimate=0)

# create a handle to watch the probability at nodes -5,0,1:
walk.watch([0,1,-5])
# create a handler to watch the entanglement
walk.watch(None,watchtype='entanglement',verbose=False)

# Propagate the CTQW using the Chebyshev method
# for t=100s in timesteps of dt=0.1
for i in np.arange(0.1,t+0.1,0.1):
	walk.propagate(i,method='chebyshev')

# plot the marginal probabilities
# after propagation over all nodes
walk.plot('out/2p_line_plot.png')

# plot the probability over time for the watched nodes
walk.plotNodes('out/2p_line_nodes.png')

# plot the entanglement over time
walk.plotEntanglement('out/2p_line_ent.png')

# export the final state
walk.exportState('out/2p_line_state.bin','bin')

# destroy the quantum walk
walk.destroy()

3P_line.py

Description

This example propagates a 3 particle continuous-time quantum walk on an infinite line

Amongst the features used, it illustrates:
  • recieving command line options using PETSc

  • the use of the chebyshev algorithm
    • setting the EigSolver tolerance, as well as the minimum eigenvalue
  • adding a diagonal defects to various nodes

  • same-node interactions between particles

  • probability vs node plot

  • Exporting the final state to a PETSc binary vector file

Output
_images/3p_line_plot.png
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 20)
t = OptDB.getReal('t', 2)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '3P Line\n'

# initialise an N (default 20) node graph CTQW
walk = qw.Line3P(N)

# Create a Hamiltonian with 2P very strong interaction.
walk.createH(interaction=10.)

# create the initial state |0,4,4>
init_state = [[0,4,4,1.0]]
walk.createInitState(init_state)

# set the eigensolver properties.
walk.EigSolver.setEigSolver(tol=1.e-2)
# underestimate the minimum eigenvalue
walk.EigSolver.setEigSolver(emin_estimate=0)

# Propagate the CTQW using the Chebyshev method for t=2s
walk.propagate(t,method='chebyshev')

# plot the marginal probabilities
# after propagation over all nodes
walk.plot('out/3p_line_plot.png')

# destroy the quantum walk
walk.destroy()

1P_3cayley.py

Description
3-cayley-graph

This example propagates a 1 particle continuous-time quantum walk on a 3-cayley tree.

Amongst the features used, it illustrates:
  • the use of the chebyshev algorithm, with eigenvalues pre-set

  • adding a diagonal defect of amplitude 3 to node 0.This acts to alter the original Hamiltonian, by the addition of a diagonal ‘disorder’ matrix \(\Gamma\); \(H = H_0 + \Gamma\) where \(\Gamma = 3|0\rangle\langle 0|\).

  • creating node handles to watch the probability at specified nodes

  • various plotting abilities:
    • probability vs node plots
    • probability vs time plots
    • graph plots
Output
_images/1p_3cayley_graph.png _images/1p_3cayley_nodes.png _images/1p_3cayley_plot.png
Required Files
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 10)
t = OptDB.getReal('t', 20)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '1P 3-cayley Tree CTQW\n'

# initialise a 10 node graph CTQW
walk = qw.Graph(N)

# Create a Hamiltonian from a file.
# A defect of amplitude 3 has been added to node 0
walk.createH('../graphs/cayley/3-cayley.txt','txt',d=[0],amp=[3.])

# create the initial state (1/sqrt(2)) (|0>+i|1>)
init_state = [[0,1.0/np.sqrt(2.0)], [1,1.0j/np.sqrt(2.0)]]
walk.createInitState(init_state)

# set the eigensolver properties. Note that
# Emin is exact, whereas and Emax has been over-estimated.
walk.EigSolver.setEigSolver(tol=1.e-2,verbose=False,emin_estimate=0.,emax_estimate=10.)

# create a handle to watch the probability at nodes 0-4:
walk.watch([0,1,2,3,4])

# Propagate the CTQW using the Chebyshev method
# for t=20s in timesteps of dt=0.01
for dt in np.arange(0.01,t+0.01,0.01):
	walk.propagate(dt,method='chebyshev')

# plot the marginal probabilities
# after propagation over all nodes
walk.plot('out/1p_3cayley_plot.png')

# plot a 3D graph representation with bars
# representing the marginal probabilities
walk.plotGraph(output='out/1p_3cayley_graph.png')

# plot the probability over time for the watched nodes
walk.plotNodes('out/1p_3cayley_nodes.png')

# destroy the quantum walk
walk.destroy()

2P_3cayley.py

Description
3-cayley-graph

This example propagates a 2 particle continuous-time quantum walk on a 3-cayley tree.

Amongst the features used, it illustrates:
  • recieving command line options using PETSc

  • the use of the chebyshev algorithm, with minimum eigenvalue pre-set

  • adding a diagonal defect to various nodes

  • same-node interactions between particles

  • creating node handles to watch the probability at specified nodes

  • creating entanglement handles to watch the entanglement

  • various plotting abilities:
    • probability vs node plots
    • probability vs time plots
    • graph plots
    • entanglement plots
  • outputs the partial trace of the density matrix in binary form

Output
_images/2p_3cayley_graph.png _images/2p_3cayley_plot.png _images/2p_3cayley_nodes_particle1.png _images/2p_3cayley_nodes_particle2.png _images/2p_3cayley_node1.png _images/2p_3cayley_ent.png
Required Files
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 10)
t = OptDB.getReal('t', 5)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '2P 3-cayley Tree CTQW\n'

# initialise a 10 node graph CTQW
walk = qw.Graph2P(N)

# Create a 2 particle interacting Hamiltonian from a file.
# A defect of amplitude 2 and 0.5 has
# been added to nodes 3 and 4 respectively.
d = [3,4]
amp = [2.0,1.5]
walk.createH('../graphs/cayley/3-cayley.txt','txt',d=d,amp=amp,layout='spring',interaction=0.5)

# create the initial state (1/sqrt(2)) (|0,1>+i|1,1>)
init_state = [[0,1,1.0/np.sqrt(2.0)], [1,1,1.0j/np.sqrt(2.0)]]
walk.createInitState(init_state)

# set the eigensolver properties.
#Note that Emin has been set exactly.
walk.EigSolver.setEigSolver(tol=1.e-2,verbose=False,emin_estimate=0.)

# create a handle to watch the probability at nodes 0-4,9:
walk.watch([0,1,2,3,4,9])
# create a handle to watch the entanglement
walk.watch(None,watchtype='entanglement',verbose=False,esolver='lapack',max_it=100)

# Propagate the CTQW using the Chebyshev method
# for t=5s in timesteps of dt=0.01
for dt in np.arange(0.01,t+0.01,0.01):
	walk.propagate(dt,method='chebyshev')

# plot the p1 and p2 marginal probabilities
# after propagation over all nodes
walk.plot('out/2p_3cayley_plot.png')

# plot a 3D graph representation with bars
# representing the p1 and p2 marginal probabilities
walk.plotGraph(output='out/2p_3cayley_graph.png')

# plot the p1 and p2 probability over time for watched node1
walk.plotNode('out/2p_3cayley_node1.png',1)
# plot the particle 1 probability over all watched nodes
walk.plotNodes('out/2p_3cayley_nodes_particle1.png',p=1)
# plot the particle 2 probability over all watched nodes
walk.plotNodes('out/2p_3cayley_nodes_particle2.png',p=2)

# plot the entanglement vs. time
walk.plotEntanglement('out/2p_3cayley_ent.png')

# export the partial trace
walk.exportPartialTrace('out/3-cayley-2p-rhoX.txt','txt',p=1)

# destroy the quantum walk
walk.destroy()

3P_3cayley.py

Description
_images/3cayley2.png

This example propagates a 3 particle continuous-time quantum walk on a 3-cayley tree.

Amongst the features used, it illustrates:
  • recieving command line options using PETSc

  • the use of the chebyshev algorithm, with minimum eigenvalue pre-set

  • adding a diagonal defect to various nodes

  • same-node interactions between particles

  • creating node handles to watch the probability at specified nodes

  • creating entanglement handles to watch the entanglement

  • various plotting abilities:
    • probability vs node plots
    • probability vs time plots
    • graph plots
    • entanglement plots
Output
_images/3p_3cayley_plot.png _images/3p_3cayley_nodes_particle1.png _images/3p_3cayley_nodes_particle3.png _images/3p_3cayley_node1.png
Required Files
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# enable command line arguments -t and -N
OptDB = PETSc.Options()
N = OptDB.getInt('N', 10)
t = OptDB.getReal('t', 5)

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '3P 3-cayley Tree CTQW\n'

# initialise a 10 node graph CTQW
walk = qw.Graph3P(N)

# Create a 3 particle interacting Hamiltonian from a file.
# A defect of amplitude 2 and 0.5 has
# been added to nodes 3 and 4 respectively.
d = [3,4]
amp = [2.0,1.5]
walk.createH('../graphs/cayley/3-cayley.txt','txt',d=d,amp=amp,layout='spring',interaction=1.5)

# create the initial state (1/sqrt(2)) (|0,1,4>+i|1,1,4>)
init_state = [[0,1,4,1.0/np.sqrt(2.0)], [1,1,4,1.j/np.sqrt(2.0)]]
walk.createInitState(init_state)

# set the eigensolver properties.
#Note that Emin has been set exactly.
walk.EigSolver.setEigSolver(tol=1.e-2,verbose=False,emin_estimate=0.)

# create a handle to watch the probability at nodes 0-4,9:
walk.watch([0,1,2,3,4,9])

# Propagate the CTQW using the Chebyshev method
# for t=5s in timesteps of dt=0.01
for dt in np.arange(0.01,t+0.01,0.01):
	walk.propagate(dt,method='chebyshev')

# plot the p1, p2, p3 marginal probabilities
# after propagation over all nodes
walk.plot('out/3p_3cayley_plot.png')

# plot the p1, p2, p3 probability over time for watched node1
walk.plotNode('out/3p_3cayley_node1.png',1)
# plot the particle 1 probability over all watched nodes
walk.plotNodes('out/3p_3cayley_nodes_particle1.png',p=1)
# plot the particle 2 probability over all watched nodes
walk.plotNodes('out/2p_3cayley_nodes_particle2.png',p=2)
# plot the particle 3 probability over all watched nodes
walk.plotNodes('out/3p_3cayley_nodes_particle3.png',p=3)

# destroy the quantum walk
walk.destroy()

1P_sr.py

Description

This example propagates a 1 particle continuous-time quantum walk on the first member of the strongly regular \((25,12,5,6)\) family of graphs.

Amongst the features used, it illustrates:
  • the use of the Krylov propagation algorithm

  • various plotting abilities:
    • probability vs node plot
    • probability and graph plot
Output
_images/1p_sr_graph.png _images/1p_sr_plot.png
Required Files
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '1P Strong Regular\n'

# initialise a 25 node graph CTQW, 
# and create a Hamiltonian from a file.
walk = qw.Graph(25)
walk.createH('../graphs/strong-regular-25-12-5-6/1.txt','txt',layout='circle')

# create the initial state |0,1>
init_state = [[0,1]]
walk.createInitState(init_state)

# Propagate the CTQW using Krylov subspace methods for t=100s
walk.propagate(100,method='krylov')

# plot the marginal probabilities
# after propagation over all nodes
walk.plot('out/1p_sr_plot.png')

# plot a 3D graph representation with bars
# representing the marginal probabilities
walk.plotGraph(output='out/1p_sr_graph.png')

# destroy the quantum walk
walk.destroy()

2P_sr.py

Description

This example propagates a 2 particle continuous-time quantum walk on the first member of the strongly regular \((25,12,5,6)\) family of graphs.

Amongst the features used, it illustrates:
  • the use of the Krylov propagation algorithm

  • various plotting abilities:
    • probability vs node plot
    • probability and graph plot
Output
_images/2p_sr_graph.png _images/2p_sr_plot.png
Required Files
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

if rank == 0:
	print '2P Strong Regular CTQW'

# initialise a 25 node 2P graph CTQW, 
# and create a Hamiltonian from a file.
walk = qw.Graph2P(25)
walk.createH('../graphs/strong-regular-25-12-5-6/1.txt','txt',layout='circle')

# create the initial state (1/sqrt(2)) (|0,1>+i|1,1>)
init_state = [[0,1,1.0/np.sqrt(2.0)], [1,1,1.0j/np.sqrt(2.0)]]
walk.createInitState(init_state)

# Propagate the CTQW using Krylov subspace methods
walk.propagate(100,method='krylov')

# plot the p1 and p2 marginal probabilities
# after propagation over all nodes
walk.plot('out/2p_sr_plot.png')

# plot a 3D graph representation with bars
# representing p1 and p2 marginal probabilities
walk.plotGraph(output='out/2p_sr_graph.png')

# destroy the quantum walk
walk.destroy()

Graph Isomorphism Examples

GraphISO_3cayley.py

Description

This example script returns the GI certificates, and also checks isomorphism of various graphs in the strongly regular \((25,12,5,6)\) family.

  1. Firstly, two permutations of the 3-Cayley tree are formed by reordering the rows and columns in such a way as to preserve the symmetry of the graph. For example,

    \[\begin{align} &\text{i)}\text{ swap rows }2\text{ & }3: \left[\begin{matrix} a & b & c\\ b & d & e\\ c & e & f \end{matrix}\right] \Rightarrow \left[\begin{matrix} b & d & e\\ a & b & c\\ c & e & f \end{matrix}\right]\\[12pt] &\text{ii)} \text{ swap columns }2\text{ & }3: \left[\begin{matrix} b & d & e\\ a & b & c\\ c & e & f \end{matrix}\right] \Rightarrow \left[\begin{matrix} d & b & e\\ b & a & c\\ e & c & f \end{matrix}\right] \end{align}\]

    – for adjacency matrices, this is process is equivalent to ‘swapping’ the label of nodes 1 and 2, and thus is isomorphic to the original adjacency matrix.

    The script then uses pyCTQW.MPI.GraphISO.AllIsomorphicQ() to test the isomorphism of all combinations of these two matrices (i.e. 1 with 2, 1 with 1, 2 with 2, and 2 with 1), and constructs a matrix containing \(1\) for isomorphic graphs, and \(0\) otherwise.

    As these graph are isomorphic, the result should be a matrix where all elements are \(1\).

  2. The second section of the program uses pyCTQW.MPI.GraphISO.AllIsomorphicQ() to test the isomorphism of all combinations of a 3-Cayley graph and a variant of the 3-Cayley graph, where an additional edge is placed between nodes 1 and 9.

    As these graph are not isomorphic, the result should be a \(2\times 2\): identity matrix.

Output
1
2
3
4
5
6
7
8
$ mpirun -np 2 GraphISO_3cayley.py
1) Testing isomorphism of all permutations:
[[ 1.  1.]
 [ 1.  1.]]

2) Testing isomorphism of all permutations:
[[ 1.  0.]
 [ 0.  1.]]
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

# create a graph isomorphism object
gi = qw.GraphISO()

# create a comparison table of two 3-cayley permutations
# present in the '../graphs/cayley/3-cayley-permutations'
# folder. As they are permutations, they are isomorphic,
# so the result should be a 2x2 matrix composed of 1.
comparisonTable = gi.AllIsomorphicQ('../graphs/cayley/3-cayley-permutations',info=False)

if rank==0:
	print '1) Testing isomorphism of all pairings:'
	print comparisonTable


# create a comparison table of a 3-cayley graph, and a
# a 3-cayley graph with an additional edge added.
# These are *not* isomorphic,
# so the result should be a 2x2 matrix identity matrix.
comparisonTable = gi.AllIsomorphicQ('../graphs/cayley/3-cayley-variant',info=False)

if rank==0:
	print '\n2) Testing isomorphism of all pairings:'
	print comparisonTable

GraphISO_sr.py

Description

This example script returns the GI certificates, and also checks isomorphism of various graphs in the strongly regular \((25,12,5,6)\) family.

  1. Firstly, two permutations of graph 1 in the SR family are formed by reordering the rows and columns in such a way as to preserve the symmetry of the graph. For example,

    \[\begin{align} &\text{i)}\text{ swap rows }2\text{ & }3: \left[\begin{matrix} a & b & c\\ b & d & e\\ c & e & f \end{matrix}\right] \Rightarrow \left[\begin{matrix} b & d & e\\ a & b & c\\ c & e & f \end{matrix}\right]\\[12pt] &\text{ii)} \text{ swap columns }2\text{ & }3: \left[\begin{matrix} b & d & e\\ a & b & c\\ c & e & f \end{matrix}\right] \Rightarrow \left[\begin{matrix} d & b & e\\ b & a & c\\ e & c & f \end{matrix}\right] \end{align}\]

    – for adjacency matrices, this is process is equivalent to ‘swapping’ the label of nodes 1 and 2, and thus is isomorphic to the original adjacency matrix.

    The GI certificates of both permutations of graph 1 are then calculated and displayed side by side. For each graph, these consist of a column of probabilities resulting from performing an interacting two particle quantum walk over all possible initial bosonic edge states; i.e.

    \[|\psi(0)\rangle = \frac{1}{\sqrt{2}}\left(|j\rangle\otimes|j\rangle + |k\rangle\otimes|k\rangle \right)~~~~\forall j,k\in\{0,1,\dots,N-1\},\]

    coupled with the frequency with which they appear. The quantum walk is performed for time \(t=2N\), allowing the walk time to propagate to all nodes. Therefore, as the graphs are isomorphic, the GI certificates should be identical.

    Finally, pyCTQW.MPI.GraphISO.isomorphicQ() is used to verify that the two graphs are in fact isomorphic.

  2. The second section of the program uses pyCTQW.MPI.GraphISO.isomorphicQ() to check whether graphs 1 and 11 of the SR family are isomorphic. They are not, and this is correctly returned by the program.

Output
 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
$ mpirun -np 2 GraphISO_sr.py

The GI certificates of two permutations of SR(25,12,5,6) graph 1:
[  1.14200170e-01   1.88140000e+04] [  1.14200170e-01   1.88140000e+04]
[  1.53741042e-01   1.54380000e+04] [  1.53741042e-01   1.54380000e+04]
[  1.76558859e-01   1.02600000e+04] [  1.76558859e-01   1.02600000e+04]
[  1.28264401e-01   9.89800000e+03] [  1.28264401e-01   9.89800000e+03]
[  1.02738389e-01   9.75000000e+03] [  1.02738389e-01   9.75000000e+03]
[  2.22889534e-01   5.45200000e+03] [  2.22889534e-01   5.45200000e+03]
[  1.43460570e-01   5.32600000e+03] [  1.43460570e-01   5.32600000e+03]
[  2.08976427e-01   4.86000000e+03] [  2.08976427e-01   4.86000000e+03]
[  1.88334629e-01   4.25600000e+03] [  1.88334629e-01   4.25600000e+03]
[  9.23774555e-02   2.46800000e+03] [  9.23774555e-02   2.46800000e+03]
[  1.98534713e-01   1.96000000e+03] [  1.98534713e-01   1.96000000e+03]
[  1.64488997e-01   1.79000000e+03] [  1.64488997e-01   1.79000000e+03]
[  2.38430321e-01   1.76200000e+03] [  2.38430321e-01   1.76200000e+03]
[  2.52977006e-01   5.34000000e+02] [  2.52977006e-01   5.34000000e+02]
[  2.65853983e-01   4.04000000e+02] [  2.65853983e-01   4.04000000e+02]
[  7.43487849e-02   3.01000000e+02] [  7.43487849e-02   3.01000000e+02]
[   0.27720404  178.        ] [   0.27720404  178.        ]

Testing isomorphism of the two permutations of SR(25,12,5,6) graph 1:
True

Testing isomorphism of graphs 1 and 11 of the SR(25,12,5,6) family:
False
Source Code

[Download source code]

 1
 2
 3
 4
 5
 6
 7
 8
 9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
#!/usr/bin/env python2.7
# initialize PETSc
import sys, petsc4py
petsc4py.init(sys.argv)
from petsc4py import PETSc
import numpy as np

# import pyCTQW as qw
import pyCTQW.MPI as qw

# get the MPI rank
rank =  PETSc.Comm.Get_rank(PETSc.COMM_WORLD)

# create a graph isomorphism object
gi = qw.GraphISO()

#=========================
#  Two Isomorphic Graphs
#=========================
# import two adjacency matrices; these are permutations
# of the first graph in the SR-25-12-5-6, and thus are isomorphic
adj1 = np.genfromtxt('../graphs/strong-regular-25-12-5-6/1-permutations/SR1_perm_1.txt')
adj2 = np.genfromtxt('../graphs/strong-regular-25-12-5-6/1-permutations/SR1_perm_3.txt')

# generate certificates for these two graphs
cert1 = gi.GIcert(adj1)
cert2 = gi.GIcert(adj2)

# print the certificates side by side;
# they should be identical
if rank == 0:
	print 'The GI certificates of two permutations of SR(25,12,5,6) graph 1:'
	for i in range(len(cert1)):
		print cert1[i], cert2[i]

if rank == 0:
	print '\nTesting isomorphism of the two permutations of SR(25,12,5,6) graph 1:'
# verify using the inbuilt checking method
checkISO = gi.isomorphicQ(adj1,adj2)
if rank == 0: print checkISO

#=========================
#Two Non-Isomorphic Graphs
#=========================
if rank == 0:
	print '\nTesting isomorphism of graphs 1 and 11 of the SR(25,12,5,6) family:'

adj1 = np.genfromtxt('../graphs/strong-regular-25-12-5-6/1.txt')
adj2 = np.genfromtxt('../graphs/strong-regular-25-12-5-6/11.txt')

checkISO = gi.isomorphicQ(adj1,adj2)
if rank == 0: print checkISO

Fortran Examples

exampleMPI.F90

Description

This Fortran example propagates a 1 particle continuous-time quantum walk on an infinite line

Amongst the features used, it illustrates:
  • calling the libctqwMPI library from Fortran
  • recieving command line options using PETSc
  • creating PETSc log stages
  • the use of the chebyshev algorithm
  • adding a diagonal defects to various nodes
  • viewing vectors via PetscViewer
Compiling and executing

To compile this example, open a terminal in the top directory of the pyCTQW-X.Y folder, and run

$ make fortran
$ make examples

The example can then be run by changing into the examples/fortran folder, and running

$ mpirun -np X exampleMPI

where X is the number of MPI processes to run.

Output
$ mpirun -np 4 exampleMPI
Vector Object: 4 MPI processes
  type: mpi
Process [0]
1.30318e-15 + 9.41288e-16 i
4.26061e-15 - 5.83926e-15 i
-2.6771e-14 - 1.97231e-14 i
-8.90705e-14 + 1.19677e-13 i
5.22435e-13 + 3.92973e-13 i
1.69261e-12 - 2.2254e-12 i
-9.24365e-12 - 7.11278e-12 i
-2.91412e-11 + 3.74127e-11 i
1.47434e-10 + 1.16317e-10 i
4.51971e-10 - 5.6522e-10 i
-2.10616e-09 - 1.70822e-09 i
-6.27416e-09 + 7.62086e-09 i
2.67489e-08 + 2.23733e-08 i
7.73784e-08 - 9.09727e-08 i
-2.99429e-07 - 2.59261e-07 i
-8.40536e-07 + 9.52525e-07 i
2.92435e-06 + 2.63331e-06 i
7.96061e-06 - 8.65084e-06 i
-2.46145e-05 - 2.31843e-05 i
-6.49345e-05 + 6.72304e-05 i
0.000175879 + 0.000174556 i
0.00044937 - 0.000439577 i
-0.00104657 - 0.00110507 i
-0.00258837 + 0.0023656 i
0.00505632 + 0.00575514 i
Process [1]
0.012099 - 0.0101714 i
-0.0191447 - 0.0239345 i
-0.0442924 + 0.0334726 i
0.053856 + 0.0761079 i
0.12025 - 0.0787511 i
-0.102846 - 0.172366 i
-0.219755 + 0.116923 i
0.111189 + 0.241323 i
0.214786 - 0.0830223 i
-0.0457964 - 0.132734 i
-0.0208922 + 0.0294486 i
0.0624386 - 0.0589731 i
-0.0404472 - 0.137791 i
-0.191971 - 0.0860889 i
-0.229877 + 0.138719 i
-0.0383627 + 0.243064 i
0.0596756 + 0.217395 i
0.21466 + 0.189964 i
0.261424 + 0.00431375 i
0.21683 - 0.0655021 i
0.178224 - 0.166243 i
0.0818483 - 0.179087 i
0.0335323 - 0.183573 i
-0.01719 - 0.145368 i
-0.0340458 - 0.120739 i
Process [2]
-0.0344425 - 0.0737224 i
-0.0328421 - 0.0407838 i
0.0103582 - 0.0348961 i
-0.0119036 + 0.0257429 i
-0.00586558 - 0.0254129 i
0.010251 - 0.000336933 i
-0.0147693 - 0.0158695 i
0.0286581 - 0.0184521 i
-0.00508709 + 0.00712859 i
0.0298955 - 0.0310717 i
0.026666 + 0.0413678 i
-0.011159 - 0.00816756 i
0.0425462 + 0.0359954 i
-0.0610777 + 0.0473539 i
-0.018832 - 0.0447773 i
-0.00248435 + 0.0247191 i
-0.0610373 - 0.0547772 i
0.0902541 - 0.0779196 i
0.0756292 + 0.10105 i
-0.0915928 + 0.061552 i
-0.0438766 - 0.0715585 i
0.049748 - 0.0280784 i
0.0163906 + 0.0313776 i
-0.0181928 + 0.00882673 i
-0.0044223 - 0.00978976 i
Process [3]
0.00492519 - 0.00207491 i
0.000916545 + 0.00233018 i
-0.0010417 + 0.000382831 i
-0.000151757 - 0.000441799 i
0.000178365 - 5.72714e-05 i
2.06318e-05 + 6.875e-05 i
-2.53648e-05 + 7.11154e-06 i
-2.35018e-06 - 8.97787e-06 i
3.05473e-06 - 7.45989e-07 i
2.27794e-07 + 1.00097e-06 i
-3.16392e-07 + 6.70093e-08 i
-1.90127e-08 - 9.66128e-08 i
2.85389e-08 - 5.20876e-09 i
1.37913e-09 + 8.16537e-09 i
-2.26541e-09 + 3.53182e-10 i
-8.75365e-11 - 6.10116e-10 i
1.59661e-10 - 2.10083e-11 i
4.88366e-12 + 4.06353e-11 i
-1.00669e-11 + 1.0998e-12 i
-2.39917e-13 - 2.42956e-12 i
5.71638e-13 - 5.06798e-14 i
1.03592e-14 + 1.31214e-13 i
-2.94057e-14 + 2.04654e-15 i
-3.89668e-16 - 6.42454e-15 i
1.43612e-15 - 7.35772e-17 i
Source Code

[Download source code]

program exampleMPI
    
    use ctqwMPI
    implicit none
#include <finclude/petsc.h>
    
    ! declare variables
    PetscLogStage  :: stage
    PetscLogDouble :: te10, te21,te11, te20, ts0, ts1, tc1, tc0
    PetscMPIInt    :: rank
    PetscErrorCode :: ierr
    PetscBool      :: flag
    PetscInt       :: i, j, its, n, d(2)
    PetscScalar    :: Emin, Emax, t, init_state(2,2), amp(2)
    PetscReal      :: Emin_error, Emax_error
    Mat            :: H, H2
    Vec            :: psi0, psi, psix, psiy
    
    ! initialize SLEPc and PETSc
    call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
    call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
    
    ! get command line arguments
    call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-n",n,flag,ierr)
    CHKERRQ(ierr)
    if (flag .eqv. .false.) n = 100
    call PetscOptionsGetInt(PETSC_NULL_CHARACTER,"-t",t,flag,ierr)
    CHKERRQ(ierr)
    if (flag .eqv. .false.) t = 10.0
    
    ! create the Hamiltonian
    d = [3,4]
    amp = [2.0,1.5]
    call PetscLogStageRegister('Hamiltonian',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    call MatCreate(PETSC_COMM_WORLD,H,ierr)
    call hamiltonian_p1_line(H,d,amp,size(d),n)
    call PetscBarrier(H,ierr)
    call PetscLogStagePop(ierr)
    ! call MatView(H,PETSC_VIEWER_STDOUT_WORLD,ierr)
    
    ! Eigenvalue solver
    call PetscTime(te10,ierr)
    call PetscLogStageRegister('Emax',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    call min_max_eigs(H,rank,Emax,Emax_error,'max','krylovschur','null',35,1.d-2,0,.false.,ierr)
    call PetscLogStagePop(ierr)
    !if (rank==0 .and. ierr==0) write(*,*)'    ',PetscRealPart(Emax),'+-',Emax_error
    call PetscTime(te11,ierr)
    
    call PetscTime(te20,ierr)
    call PetscLogStageRegister('Emin',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    call min_max_eigs(H,rank,Emin,Emin_error,'min','krylovschur','null',35,1.d-2,0,.false.,ierr)
    call PetscLogStagePop(ierr)
    call PetscBarrier(H,ierr)
   ! if (rank==0 .and. ierr==0) write(*,*)'    ',PetscRealPart(Emin),'+-', Emin_error
    call PetscTime(te21,ierr)
    
    ! create vectors
    call VecCreate(PETSC_COMM_WORLD,psi0,ierr)
    call VecSetSizes(psi0,PETSC_DECIDE,n,ierr)
    call VecSetFromOptions(psi0,ierr)
    
    call VecDuplicate(psi0,psi,ierr)
    
    ! create initial state
    call PetscLogStageRegister('InitState',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    init_state(1,:) = [0.,1.0/sqrt(2.0)]
    init_state(2,:) = [1.,1.0/sqrt(2.0)]
    call p1_init(psi0,init_state,2,n)
    !call VecView(psi0,PETSC_VIEWER_STDOUT_WORLD,ierr)
    call PetscBarrier(psi0,ierr)
    call PetscLogStagePop(ierr)
    
    ! matrix exponential
    call PetscTime(ts0,ierr)
    call PetscLogStageRegister('SLEPc expm',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    call qw_krylov(H,t,psi0,psi)
    !call VecView(psi,PETSC_VIEWER_STDOUT_WORLD,ierr)
    call PetscBarrier(psi,ierr)
    call PetscLogStagePop(ierr)
    call PetscTime(ts1,ierr)
    
    ! QW chebyshev
    call PetscTime(tc0,ierr)
    call PetscLogStageRegister('Chebyshev',stage,ierr)
    call PetscLogStagePush(stage,ierr)
    call qw_cheby(psi0,psi,t,H,0.*PETSc_i,Emax,rank,n)
    call VecView(psi,PETSC_VIEWER_STDOUT_WORLD,ierr)
    call PetscBarrier(psi,ierr)
    call PetscLogStagePop(ierr)
    call PetscTime(tc1,ierr)

    ! destroy matrix/SLEPc
    call MatDestroy(H,ierr)
    call VecDestroy(psi,ierr)
    call VecDestroy(psi0,ierr)
    
    call PetscFinalize(ierr)

end program exampleMPI
  • Fortran and Python bindings (in the form of a library and module respectively)
  • Supports MPI through the use of the PETSc and SLEPc high-performance sparse matrix libraries (CUDA support planned)
  • Has in-built support for infinite-line Hamiltonians
  • Can import custom adjancency matrices
  • Supports one, two and three continuous-time quantum walkers (with or without interactions)
  • Import and export matrices/states in binary or text format
  • Python module supports plotting and visualisation using matplotlib and networkx
  • Entanglement calculations
  • Ability to place diagonal defects/barriers on graph vertices
  • MPI graph isomorphism methods, for both 2 and 3 interacting particles
Gallery
_images/2p_3cayley_nodes_particle1.png _images/1p_3cayley_graph.png _images/1p_sr_graph.png _images/2p_line_plot.png