Julia+GPU
This tutorial covers matrix operations on large, sparse matrices using the Julia language and the CUDA package for interacting with the GPU. The concrete example used is the time evolution of the Ising model.
Matrix creation
The following code creates the (Hamiltonian) matrix on N qubits:
using Pkg
Pkg.activate(".")
using LinearAlgebra
using SparseArrays
using ExponentialUtilities
using CUDA
using CUDA.CUSPARSE
CUDA.allowscalar(false)
function IsingHamiltonian_sparse(N)
id = [1.0 0; 0 1.] |> sparse
X = [0 1.; 1. 0] |> sparse
Z = [1.0 0; 0 -1.] |> sparse
global ham = spzeros(2^N, 2^N)
global q_term_ops = fill(id, N)
global q_term_ops[1] = X + Z
for i in 1:N
global ham += (pi/4)*foldl(kron, q_term_ops)
global q_term_ops = circshift(q_term_ops,1)
end
global qq_term_ops = fill(id, N)
global qq_term_ops[1] = Z
global qq_term_ops[2] = Z
for i in 1:N-1
global ham -= (3*pi/4)*foldl(kron, qq_term_ops)
global qq_term_ops = circshift(qq_term_ops,1)
end
return ham
end
The following code is used to study time evolution of a state whose initial configuration is the all up state on 2 qubits.
num_qubits = 2
up = [1. 0] |> sparse
state = vec(foldl(kron, fill(up,num_qubits)))
time = 1
ham = IsingHamiltonian_sparse(num_qubits)
hamCuda = CuSparseMatrixCSR(ham)
stateCuda = CuArray(state)
stateCuda = expv(time,-im*hamCuda,stateCuda)
println(stateCuda)
Environment setup
This script is saved in a file called ising_model.jl. Existing Julia modules are used to run the code on Hyak. To load the module and install the packages, do the following.
$ cd /gscratch/uwqcc/<netid>
$ mkdir julia_env
$ cd julia_env
$ srunĀ -A uwqcc -p gpu-l40 --time=1:00:00 --mem=10G --pty /bin/bash
$ module load elbert/julia/1.11.1/1.11.1
$ julia
julia> ] # enter the package manager by typing ']'
pkg> activate .
pkg> add 'ExponentialUtilities' 'CUDA' 'SparseArrays'
pkg> up
pkg> # hit backspace or esc to exit package manager
julia> exit() # exit julia
The ising_model.jl code should be placed within the newly-created julia_env environment.
Job submission
The following slurm submit file ising_model.slurm is used to submit the batch job.
#!/bin/bash
#SBATCH --job-name=ising_model
#SBATCH --mail-type=ALL # when to email about job status updates
#SBATCH --mail-user=<netid>@uw.edu # email for notifications about job status
#SBATCH -p gpu-l40 # partition
#SBATCH -A uwqcc # account
#SBATCH --nodes=1 # number of nodes
#SBATCH --time=1:00:00 # wall time
#SBATCH --ntasks=1 # number of threads
#SBATCH --gpus=l40:1 # number of gpus
#SBATCH --mem=40G # amount of RAM
#SBATCH --chdir=/gscratch/uwqcc/<netid>/julia_env
module load elbert/julia/1.11.1/1.11.1
julia ising_model.jl
exit $?
To submit the ising_model.slurm job, run
$ sbatch ising_model.slurm
Once the job is finished, the output should read:
ComplexF64[0.3315866240881443 + 0.7909525883567952im, -0.2641567658718962 - 0.20881547521215435im, -0.2641567658718962 - 0.20881547521215435im, -0.12110260105494038 - 0.15170060786722808im],
which is the time-evolved statevector of the two-qubit system.