Other extensions

Apart from ITensor extension, FuzzifiED also provides other extensions, viz. HDF5 extension and CUDA extension.

HDF5 extension

The HDF5 extension supports writing the types Confs, Basis, Terms, Operator, OpMat{ComplexF64} and OpMat{Float64} into HDF5 files and reading them from groups and subgroups in HDF5 format.

using HDF5 
h5open(file_name, "cw")
# include the file name as a string 
# Modes : "cw" for write and "r" for read
...
close(f)

To write, include in the middle

write(f, group_name :: String, cfs :: Confs)
write(f, group_name :: String, bs  :: Basis)
write(f, group_name :: String, tms :: Terms)
write(f, group_name :: String, op  :: Operator)
write(f, group_name :: String, mat :: OpMat{ComplexF64})
write(f, group_name :: String, mat :: OpMat{Float64})

To read, include in the middle

cfs = read(f, group_name :: String, Confs)
bs  = read(f, group_name :: String, Basis)
tms = read(f, group_name :: String, Terms)
op  = read(f, group_name :: String, Operator)
mat = read(f, group_name :: String, OpMat{ComplexF64})
mat = read(f, group_name :: String, OpMat{Float64})

CUDA extension

CUDA.CUSPARSE.CuSparseMatrixCSCMethod
CUSPARSE.CuSparseMatrixCSC(mat :: OpMat{ComplexF64})
CUSPARSE.CuSparseMatrixCSC(mat :: OpMat{Float64})

converts the OpMat objects to a CuSparseMatrixCSC object in the CUDA.CUSPARSE package.

source
FuzzifiED.GetEigensystemCudaMethod
GetEigensystemCuda(mat :: OpMat{ComplexF64}, nst :: Int64 ; initvec :: Vector{ComplexF64}, num_th :: Int64, disp_std :: Bool, kwargs...) :: Tuple{Vector{ComplexF64}, CuArray{ComplexF64, 2, CUDA.DeviceMemory}}
GetEigensystemCuda(mat :: OpMat{Float64}, nst :: Int64 ; initvec :: Vector{Float64}, num_th :: Int64, disp_std :: Bool, kwargs...) :: Tuple{Vector{Float64}, CuArray{Float64, 2, CUDA.DeviceMemory}}

This method uses Julia KrylovKit package to calculate the lowest eigenstates of sparse matrix. The sparse matrix multiplication is realised by CUDA.CUSPARSE. For an example, refer to ising_spectrum_cuda.jl.

Arguments

  • mat :: OpMat{ComplexF64} or mat :: OpMat{Float64} is the matrix ;
  • nst :: Int64 is the number of eigenstates to be calculated ;
  • tol :: Float64 is the tolerence for the KrylovKit process. The default value is 1E-8 ;
  • ncv :: Int64 is the maximum dimension of the Krylov subspace. The default value is max(2 * nst, nst + 10). If krylovdim is also given, ncv will not be used.
  • initvec is the initial vector. Facultative, a random initialisation CUDA.rand(T, mat.dimd) by default ;
  • num_th :: Int64, the number of threads. Facultative, NumThreads by default.
  • disp_std :: Bool, whether or not the log shall be displayed. Facultative, !SilentStd by default.
  • kwargs... is the options that will directly sent into eigsolve, see its documentation for detail.

Output

  • A length-nst array that has the same type as mat recording the eigenvalues, and
  • A dimd*nst matrix that has the same type as mat where every column records an eigenstate.
source