-
Notifications
You must be signed in to change notification settings - Fork 1
Commit
This commit does not belong to any branch on this repository, and may belong to a fork outside of the repository.
Got multidim polynomials working. Default constructed for bitmaps of …
…arb dimension. Refact tests
- Loading branch information
Showing
10 changed files
with
193 additions
and
123 deletions.
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,32 @@ | ||
using ITensorNumericalAnalysis | ||
|
||
using Graphs: SimpleGraph, uniform_tree | ||
using NamedGraphs: NamedGraph | ||
using ITensors: siteinds | ||
using Random: Random | ||
|
||
L = 12 | ||
Random.seed!(1234) | ||
g = NamedGraph(SimpleGraph(uniform_tree(L))) | ||
bit_map = BitMap(g; map_dimension=3) | ||
s = siteinds(g, bit_map) | ||
|
||
println( | ||
"Constructing the 3D function f(x,y,z) = x³(y + y²) + cosh(πz) as a tensor network on a randomly chosen tree with $L vertices", | ||
) | ||
ψ_fx = poly_itn(s, bit_map, [0.0, 0.0, 0.0, 1.0]; dimension=1) | ||
ψ_fy = poly_itn(s, bit_map, [0.0, 1.0, 1.0, 0.0]; dimension=2) | ||
ψ_fz = cosh_itn(s, bit_map; k=Float64(pi), dimension=3) | ||
ψxyz = ψ_fx * ψ_fy + ψ_fz | ||
|
||
ψxyz = truncate(ψxyz; cutoff=1e-12) | ||
println("Maximum bond dimension of the network is $(maxlinkdim(ψxyz))") | ||
|
||
x, y, z = 0.125, 0.625, 0.5 | ||
fxyz_xyz = calculate_fxyz(ψxyz, [x, y, z]) | ||
println( | ||
"Tensor network evaluates the function as $fxyz_xyz at the co-ordinate: (x,y,z) = ($x, $y, $z)", | ||
) | ||
println( | ||
"Actual value of the function is $(x^3 * (y + y^2) + cosh(pi * z)) at the co-ordinate: (x,y,z) = ($x, $y, $z)", | ||
) |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,78 @@ | ||
using ITensors: Index, ITensor, dim | ||
|
||
"""Exponent on x_i for the tensor Q(x_i) on the tree""" | ||
function f_alpha_beta(α::Vector{Int64}, beta::Int64) | ||
return !isempty(α) ? max(0, beta - sum(α)) : max(0, beta) | ||
end | ||
|
||
"""Coefficient on x_i for the tensor Q(x_i) on the tree""" | ||
function _coeff(N::Int64, α::Vector{Int64}, beta) | ||
@assert length(α) == N - 1 | ||
return if N == 1 | ||
1 | ||
else | ||
prod([binomial(f_alpha_beta(α[1:(N - 1 - i)], beta), α[N - i]) for i in 1:(N - 1)]) | ||
end | ||
end | ||
|
||
"""Constructor for the tensor that sits on a vertex of degree N""" | ||
function Q_N_tensor( | ||
N::Int64, siteind::Index, αind::Vector{Index}, betaind::Index, xivals::Vector{Float64} | ||
) | ||
@assert length(αind) == N - 1 | ||
@assert length(xivals) == dim(siteind) | ||
n = dim(betaind) - 1 | ||
@assert all(x -> x == n + 1, dim.(αind)) | ||
|
||
link_dims = [n + 1 for i in 1:N] | ||
dims = vcat([dim(siteind)], link_dims) | ||
Q_N_array = zeros(Tuple(dims)) | ||
for (i, xi) in enumerate(xivals) | ||
for j in 0:((n + 1)^(N) - 1) | ||
is = digits(j; base=n + 1, pad=N) | ||
f = f_alpha_beta(is[1:(N - 1)], last(is)) | ||
Q_N_array[(i, Tuple(is + ones(Int64, (N)))...)...] = | ||
_coeff(N, is[1:(N - 1)], last(is)) * (xi^f) | ||
end | ||
end | ||
|
||
return ITensor(Q_N_array, siteind, αind, betaind) | ||
end | ||
|
||
"""Tensor for transferring the alpha index (beta ind here) of a QN tensor (defined above) across multiple inds (alpha_inds)""" | ||
#Needed for building multi-dimensional polynomials | ||
function transfer_tensor(phys_ind::Index, beta_ind::Index, alpha_inds::Vector) | ||
virt_inds = vcat(beta_ind, alpha_inds) | ||
inds = vcat(phys_ind, virt_inds) | ||
virt_dims = dim.(virt_inds) | ||
@assert allequal(virt_dims) | ||
dims = vcat(dim(phys_ind), virt_dims) | ||
N = length(alpha_inds) | ||
T_array = zeros(Tuple(dims)) | ||
for i in 0:(dim(phys_ind) - 1) | ||
for j in 0:(dim(beta_ind) - 1) | ||
if !isempty(alpha_inds) | ||
for k in 0:((first(virt_dims))^(N) - 1) | ||
is = digits(k; base=first(virt_dims), pad=N) | ||
if sum(is) == j | ||
T_array[(i + 1, j + 1, Tuple(is + ones(Int64, (N)))...)...] = 1.0 | ||
end | ||
end | ||
else | ||
T_array[i + 1, j + 1] = j == 0 ? 1 : 0 | ||
end | ||
end | ||
end | ||
return ITensor(T_array, phys_ind, beta_ind, alpha_inds) | ||
end | ||
|
||
"""Given a tree find the edge coming from the vertex v which is directed towards `source_vertex`""" | ||
function get_edge_toward_vertex(g::AbstractGraph, v, source_vertex) | ||
for vn in neighbors(g, v) | ||
if length(a_star(g, vn, source_vertex)) < length(a_star(g, v, source_vertex)) | ||
return NamedEdge(v => vn) | ||
end | ||
end | ||
|
||
return error("Couldn't find edge. Perhaps graph is not a tree or not connected.") | ||
end |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,11 @@ | ||
@eval module $(gensym()) | ||
using ITensorNumericalAnalysis: ITensorNumericalAnalysis | ||
using Test: @testset | ||
|
||
@testset "Test examples" begin | ||
example_files = ["2d_laplace_solver.jl", "construct_multi_dimensional_function.jl"] | ||
@testset "Test $example_file" for example_file in example_files | ||
include(joinpath(pkgdir(ITensorNumericalAnalysis), "examples", example_file)) | ||
end | ||
end | ||
end |
Oops, something went wrong.