Skip to content

Commit be1e2a1

Browse files
leburgellkdvos
andauthored
Expand exact_diagonalization and FiniteMPS docstrings (#261)
* Expand `exact_diagonalization` and `FiniteMPS` docstrings * Language * Update src/algorithms/ED.jl Co-authored-by: Lukas Devos <[email protected]> * And some more `opp -> H` --------- Co-authored-by: Lukas Devos <[email protected]>
1 parent c322570 commit be1e2a1

File tree

2 files changed

+40
-17
lines changed

2 files changed

+40
-17
lines changed

src/algorithms/ED.jl

Lines changed: 36 additions & 15 deletions
Original file line numberDiff line numberDiff line change
@@ -1,24 +1,45 @@
11
"""
2-
exact_diagonalization(opp::FiniteMPOHamiltonian;
3-
sector=first(sectors(oneunit(physicalspace(opp, 1)))),
4-
len::Int=length(opp), num::Int=1, which::Symbol=:SR,
2+
exact_diagonalization(H::FiniteMPOHamiltonian;
3+
sector=first(sectors(oneunit(physicalspace(H, 1)))),
4+
len::Int=length(H), num::Int=1, which::Symbol=:SR,
55
alg=Defaults.alg_eigsolve(; dynamic_tols=false))
66
-> vals, state_vecs, convhist
77
88
Use [`KrylovKit.eigsolve`](@extref) to perform exact diagonalization on a
99
`FiniteMPOHamiltonian` to find its eigenvectors as `FiniteMPS` of maximal rank, essentially
1010
equivalent to dense eigenvectors.
11+
12+
### Arguments
13+
- `H::FiniteMPOHamiltonian`: the Hamiltonian to diagonalize.
14+
15+
### Keyword arguments
16+
- `sector=first(sectors(oneunit(physicalspace(H, 1))))`: the total charge of the
17+
eigenvectors, which is chosen trivial by default.
18+
- `len::Int=length(H)`: the length of the system.
19+
- `num::Int=1`: the number of eigenvectors to find.
20+
- `which::Symbol=:SR`: the kind eigenvalues to find, see [`KrylovKit.eigsolve`](@extref).
21+
- `alg=Defaults.alg_eigsolve(; dynamic_tols=false)`: the diagonalization algorithm to use,
22+
see [`KrylovKit.eigsolve`](@extref).
23+
24+
!!! note "Valid `sector` values"
25+
The total charge of the eigenvectors is imposed by adding a charged auxiliary space as
26+
the leftmost virtualspace of each eigenvector. Specifically, this is achieved by passing
27+
`left=Vect[typeof(sector)](sector => 1)` to the [`FiniteMPS`](@ref) constructor. As
28+
such, the only valid `sector` values (i.e. `sector` values for which the corresponding
29+
eigenstates have valid fusion channels) are those that occur in the dual of the fusion
30+
of all the physical spaces in the system.
31+
1132
"""
12-
function exact_diagonalization(opp::FiniteMPOHamiltonian;
13-
sector=first(sectors(oneunit(physicalspace(opp, 1)))),
14-
len::Int=length(opp), num::Int=1, which::Symbol=:SR,
33+
function exact_diagonalization(H::FiniteMPOHamiltonian;
34+
sector=first(sectors(oneunit(physicalspace(H, 1)))),
35+
len::Int=length(H), num::Int=1, which::Symbol=:SR,
1536
alg=Defaults.alg_eigsolve(; dynamic_tols=false))
1637
left = ℂ[typeof(sector)](sector => 1)
1738
right = oneunit(left)
1839

1940
middle_site = Int(round(len / 2))
2041

21-
Ot = eltype(opp[1])
42+
Ot = eltype(H[1])
2243

2344
mpst_type = tensormaptype(spacetype(Ot), 2, 1, storagetype(Ot))
2445
mpsb_type = tensormaptype(spacetype(Ot), 1, 1, storagetype(Ot))
@@ -28,27 +49,27 @@ function exact_diagonalization(opp::FiniteMPOHamiltonian;
2849
ACs = Vector{Union{Missing,mpst_type}}(missing, len)
2950

3051
for i in 1:(middle_site - 1)
31-
ALs[i] = isomorphism(storagetype(Ot), left * physicalspace(opp, i),
32-
fuse(left * physicalspace(opp, i)))
52+
ALs[i] = isomorphism(storagetype(Ot), left * physicalspace(H, i),
53+
fuse(left * physicalspace(H, i)))
3354
left = _lastspace(ALs[i])'
3455
end
3556
for i in len:-1:(middle_site + 1)
3657
ARs[i] = _transpose_front(isomorphism(storagetype(Ot),
37-
fuse(right * physicalspace(opp, i)'),
38-
right * physicalspace(opp, i)'))
58+
fuse(right * physicalspace(H, i)'),
59+
right * physicalspace(H, i)'))
3960
right = _firstspace(ARs[i])
4061
end
41-
ACs[middle_site] = randomize!(similar(opp[1][1, 1, 1, 1],
42-
left * physicalspace(opp, middle_site) right))
62+
ACs[middle_site] = randomize!(similar(H[1][1, 1, 1, 1],
63+
left * physicalspace(H, middle_site) right))
4364
norm(ACs[middle_site]) == 0 && throw(ArgumentError("invalid sector"))
4465
normalize!(ACs[middle_site])
4566

4667
#construct the largest possible finite mps of that length
4768
state = FiniteMPS(ALs, ARs, ACs, Cs)
48-
envs = environments(state, opp)
69+
envs = environments(state, H)
4970

5071
#optimize the middle site. Because there is no truncation, this single site captures the entire possible hilbert space
51-
H_ac = ∂∂AC(middle_site, state, opp, envs) # this linear operator is now the actual full hamiltonian!
72+
H_ac = ∂∂AC(middle_site, state, H, envs) # this linear operator is now the actual full hamiltonian!
5273
(vals, vecs, convhist) = eigsolve(H_ac, state.AC[middle_site], num, which, alg)
5374

5475
state_vecs = map(vecs) do v

src/states/finitemps.jl

Lines changed: 4 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -34,7 +34,9 @@ By convention, we have that:
3434
FiniteMPS(As::Vector{<:GenericMPSTensor}; normalize=false, overwrite=false)
3535
3636
Construct an MPS via a specification of physical and virtual spaces, or from a list of
37-
tensors `As`. All cases reduce to the latter.
37+
tensors `As`. All cases reduce to the latter. In particular, a state with a non-trivial
38+
total charge can be constructed by passing a non-trivially charged vector space as the
39+
`left` or `right` virtual spaces.
3840
3941
### Arguments
4042
- `As::Vector{<:GenericMPSTensor}`: vector of site tensors
@@ -50,7 +52,7 @@ tensors `As`. All cases reduce to the latter.
5052
- `maxvirtualspace::S`: maximum virtual space
5153
5254
### Keywords
53-
- `normalize`: normalize the constructed state
55+
- `normalize=true`: normalize the constructed state
5456
- `overwrite=false`: overwrite the given input tensors
5557
- `left=oneunit(S)`: left-most virtual space
5658
- `right=oneunit(S)`: right-most virtual space

0 commit comments

Comments
 (0)