Other extensions
Apart from ITensor extension, FuzzifiED also provides other extensions, viz. HDF5 extension, KrylovKit extension, CUDA extension and SparseArrays 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. This extension requires the packages HDF5
. To use this extension, include at the heading
using HDF5
A typical file operation process looks like
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})
SparseArrays extension
The SparseArrays extension supports the conversion between OpMat
in FuzzifiED and the SparseMatrixCSC
and Matrix
format. This extension requires the packages SparseArrays
. To use this extension, include at the heading
using SparseArrays
SparseArrays.SparseMatrixCSC
— MethodSparseMatrixCSC(mat :: OpMat{ComplexF64}) :: SparseMatrixCSC{Int64,ComplexF64}
SparseMatrixCSC(mat :: OpMat{Float64}) :: SparseMatrixCSC{Int64,Float64}
converts the OpMat
objects to a SparseMatrixCSC
object in the SparseArrays
package.
Base.Matrix
— MethodMatrix(mat :: OpMat{ComplexF64}) :: Matrix{ComplexF64}
Matrix(mat :: OpMat{Float64}) :: Matrix{Float64}
converts the OpMat
objects to a full matrix.
FuzzifiED.OpMat
— MethodOpMat(matcsc :: SparseMatrixCSC{Int64,ComplexF64}) :: OpMat{ComplexF64}
OpMat(matcsc :: SparseMatrixCSC{Int64,Float64}) :: OpMat{Float64}
converts the SparseMatrixCSC
object in the SparseArrays
package to an OpMat
objects.
KrylovKit extension
The KrylovKit extension supports an interface diagonalising the sparse matrix using the KrylovKit
pakage in Julia. This extension requires the packages KrylovKit
. To use this extension, include at the heading
using KrylovKit
Besides Arpack in Fortran, we also provide an interface calling KrylovKit in Julia.
FuzzifiED.GetEigensystemKrylov
— MethodGetEigensystemKrylov(mat :: OpMat{ComplexF64}, nst :: Int64 ; initvec :: Vector{ComplexF64}, num_th :: Int64, disp_std :: Bool, kwargs...) :: Tuple{Vector{ComplexF64}, Matrix{ComplexF64}}
GetEigensystemKrylov(mat :: OpMat{Float64}, nst :: Int64 ; initvec :: Vector{Float64}, num_th :: Int64, disp_std :: Bool, kwargs...) :: Tuple{Vector{Float64}, Matrix{Float64}}
This method calls the eigsolve
from Julia KrylovKit.jl
package instead of Arpack from Fortran to calculate the lowest eigenstates of sparse matrix. The performance should be similar. For an example, refer to ising_spectrum_krylov.jl
.
Arguments
mat :: OpMat{ComplexF64}
ormat :: 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 is1E-8
.ncv :: Int64
is the maximum dimension of the Krylov subspace. The default value ismax(2 * nst, nst + 10)
. Ifkrylovdim
is also given,ncv
will not be used.initvec :: Vector{ComplexF64}
orinitvec :: Vector{Float64}
is the initial vector. Facultative, a random initialisation 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 intoeigsolve
, see its documentation for detail.
Output
- A length-
nst
array that has the same type asmat
recording the eigenvalues, and - A
dimd
*nst
matrix that has the same type asmat
where every column records an eigenstate.
CUDA extension
The CUDA extension supports the conversion between OpMat
in FuzzifiED and the CuSparseMatrixCSC
in CUDA.CUSPARSE
as well as the acceleration of diagonalisation on GPU. This extension requires the packages CUDA
, KrylovKit
and SparseArrays
. To use this extension, include at the heading
using CUDA, KrylovKit, SparseArrays
CUDA.CUSPARSE.CuSparseMatrixCSC
— MethodCUSPARSE.CuSparseMatrixCSC(mat :: OpMat{ComplexF64})
CUSPARSE.CuSparseMatrixCSC(mat :: OpMat{Float64})
converts the OpMat
objects to a CuSparseMatrixCSC
object in the CUDA.CUSPARSE
package.
FuzzifiED.GetEigensystemCuda
— MethodGetEigensystemCuda(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}
ormat :: 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 is1E-8
.ncv :: Int64
is the maximum dimension of the Krylov subspace. The default value ismax(2 * nst, nst + 10)
. Ifkrylovdim
is also given,ncv
will not be used.initvec
is the initial vector. Facultative, a random initialisationCUDA.rand(T, mat.dimd)
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 intoeigsolve
, see its documentation for detail.
Output
- A length-
nst
array that has the same type asmat
recording the eigenvalues, and - A
dimd
*nst
matrix that has the same type asmat
where every column records an eigenstate.