@@ -4,17 +4,22 @@ using ITensors:
44 ITensor,
55 Index,
66 QN,
7+ apply,
78 dag,
89 delta,
910 inds,
11+ mapprime,
1012 noprime,
13+ norm,
1114 op,
15+ permute,
1216 prime,
1317 random_itensor,
1418 replaceind,
1519 replaceinds,
16- sim
17- using ITensorNetworks. ITensorsExtensions: map_eigvals
20+ sim,
21+ swapprime
22+ using ITensorNetworks. ITensorsExtensions: eigendecomp, map_eigvals
1823using StableRNGs: StableRNG
1924using Test: @test , @testset
2025@testset " ITensorsExtensions" begin
@@ -55,9 +60,9 @@ using Test: @test, @testset
5560 rng = StableRNG (1234 )
5661 A = random_itensor (rng, elt, i, j)
5762 P = A * prime (dag (A), i)
58- sqrtP = map_eigvals (sqrt, P, i, i ' ; ishermitian= true )
59- inv_P = dag (map_eigvals (inv, P, i, i ' ; ishermitian= true ))
60- inv_sqrtP = dag (map_eigvals (inv ∘ sqrt, P, i, i ' ; ishermitian= true ))
63+ sqrtP = map_eigvals (sqrt, P, [i], [i ' ] ; ishermitian= true )
64+ inv_P = dag (map_eigvals (inv, P, [i], [i ' ] ; ishermitian= true ))
65+ inv_sqrtP = dag (map_eigvals (inv ∘ sqrt, P, [i], [i ' ] ; ishermitian= true ))
6166
6267 new_ind = noprime (sim (i' ))
6368 sqrtPdag = replaceind (dag (sqrtP), i' , new_ind)
@@ -72,5 +77,66 @@ using Test: @test, @testset
7277 I = replaceind (inv_sqrtP * sqrtP, new_ind, i)
7378 @test I ≈ op (" I" , i)
7479 end
80+
81+ @testset " Fermionic eigendecomp" begin
82+ ITensors. enable_auto_fermion ()
83+ s1 = Index ([QN (" Nf" , 0 , - 1 )=> 2 , QN (" Nf" , 1 , - 1 )=> 2 ], " Site,Fermion,n=1" )
84+ s2 = Index ([QN (" Nf" , 0 , - 1 )=> 2 , QN (" Nf" , 1 , - 1 )=> 2 ], " Site,Fermion,n=2" )
85+
86+ # Make a random Hermitian matrix-like 4th order ITensor
87+ T = random_itensor (s1' , s2' , dag (s2), dag (s1))
88+ T = apply (T, swapprime (dag (T), 0 => 1 ))
89+ @test T ≈ swapprime (dag (T), 0 => 1 ) # check Hermitian
90+
91+ Ul, D, Ur = eigendecomp (T, [s1' , s2' ], [dag (s1), dag (s2)]; ishermitian= true )
92+
93+ @test Ul* D* Ur ≈ T
94+ ITensors. disable_auto_fermion ()
95+ end
96+
97+ @testset " Fermionic map eigvals tests" begin
98+ ITensors. enable_auto_fermion ()
99+ s1 = Index ([QN (" Nf" , 0 , - 1 )=> 2 , QN (" Nf" , 1 , - 1 )=> 2 ], " Site,Fermion,n=1" )
100+ s2 = Index ([QN (" Nf" , 0 , - 1 )=> 2 , QN (" Nf" , 1 , - 1 )=> 2 ], " Site,Fermion,n=2" )
101+
102+ # Make a random Hermitian matrix ITensor
103+ M = random_itensor (s1' , dag (s1))
104+ # M = mapprime(prime(M)*swapprime(dag(M),0=>1),2=>1)
105+ M = apply (M, swapprime (dag (M), 0 => 1 ))
106+
107+ # Make a random Hermitian matrix-like 4th order ITensor
108+ T = random_itensor (s1' , s2' , dag (s2), dag (s1))
109+ T = apply (T, swapprime (dag (T), 0 => 1 ))
110+
111+ # Matrix test
112+ sqrtM = map_eigvals (sqrt, M, [s1' ], [dag (s1)]; ishermitian= true )
113+ @test M ≈ apply (sqrtM, sqrtM)
114+
115+ # # Tensor test
116+ sqrtT = map_eigvals (sqrt, T, [s1' , s2' ], [dag (s1), dag (s2)]; ishermitian= true )
117+ @test T ≈ apply (sqrtT, sqrtT)
118+
119+ # Permute and test again
120+ T = permute (T, dag (s2), s2' , dag (s1), s1' )
121+ sqrtT = map_eigvals (sqrt, T, [s1' , s2' ], [dag (s1), dag (s2)]; ishermitian= true )
122+ @test T ≈ apply (sqrtT, sqrtT)
123+
124+ # # Explicitly passing indices in different, valid orders
125+ sqrtT = map_eigvals (sqrt, T, [s2' , s1' ], [dag (s2), dag (s1)]; ishermitian= true )
126+ @test T ≈ apply (sqrtT, sqrtT)
127+ sqrtT = map_eigvals (sqrt, T, [dag (s2), dag (s1)], [s2' , s1' ], ; ishermitian= true )
128+ @test T ≈ apply (sqrtT, sqrtT)
129+ sqrtT = map_eigvals (sqrt, T, [dag (s1), dag (s2)], [s1' , s2' ], ; ishermitian= true )
130+ @test T ≈ apply (sqrtT, sqrtT)
131+
132+ # Test bosonic index case while fermion system is enabled
133+ b = Index ([QN (" Nb" , 0 )=> 2 , QN (" Nb" , 1 )=> 2 ])
134+ T = random_itensor (b' , dag (b))
135+ T = apply (T, swapprime (dag (T), 0 => 1 ))
136+ sqrtT = map_eigvals (sqrt, T, [b' ], [dag (b)]; ishermitian= true )
137+ @test T ≈ apply (sqrtT, sqrtT)
138+
139+ ITensors. disable_auto_fermion ()
140+ end
75141end
76142end
0 commit comments