/
plane_wave_excitation.py
1393 lines (1193 loc) · 62.7 KB
/
plane_wave_excitation.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""Plane wave excitations ansatz
The quasiparticle ansatz, utilizing plane wave states, offers a powerful approach for computing
low-energy excitations. For finite systems, we could employ DMRG with different charge sectors to
find these states. For infinite systems, an efficient algorithm was introduced
in :cite:`haegeman2012`. By working directly in the tangent space of a uniform
MPS :cite:`vanderstraeten2019`, we can use translational invariance to specify momentum.
The plane wave excitation finds excitation of a :class:`~tenpy.networks.uniform_mps.UniformMPS`.
By summing over all states, where one tensor is replaced with an excited tensor B, we get a state
with a fixed momentum :class:`~tenpy.networks.momentum_mps.MomentumMPS`. The tensors B are
decomposed into -B- = -V-X-, where V is the orthogonal complement of the usual A-tensors of the
uniform MPS, and X contains the variational parameters. The algorithm constructs an effective
Hamiltonian for X and finds the low-energy states using an iterative eigensolver. Since we don't
sweep as e.g. in DMRG, more Lanczos steps may be required. Increase those until the energy is
converged!
Additionally, we can specify X to be in a given charge sector.
The :class:`PlaneWaveExcitationEngine` optimizes the X for each tensor in the unit cell. The ansatz
can be extended to include excitations that span several sites. This is implemented
in :class:`MultiSitePlaneWaveExcitationEngine`. Note that with the current implementation, the
numerical costs scale exponentially with the number of exciting sites.
"""
# Copyright (C) TeNPy Developers, GNU GPLv3
import numpy as np
import logging
logger = logging.getLogger(__name__)
from ..linalg import np_conserved as npc
from ..linalg.charges import LegPipe
from ..networks.momentum_mps import MomentumMPS
from ..networks.mpo import MPOEnvironment, MPOTransferMatrix
from ..linalg.krylov_based import GMRES, LanczosGroundState, Arnoldi
from ..linalg.sparse import NpcLinearOperator, SumNpcLinearOperator, BoostNpcLinearOperator
from ..algorithms.algorithm import Algorithm
from ..algorithms.mps_common import ZeroSiteH
__all__ = [
'append_right_env', 'append_left_env', 'construct_orthogonal', 'PlaneWaveExcitationEngine',
'MultiSitePlaneWaveExcitationEngine'
]
def append_right_env(As, Bs, R, Ws=None):
"""Contract all tensors in As and Bs to the right environment R.
If Ws given, contract as a MPO environment.
Parameters
----------
As : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the top leg of the environment, labels are ``vL, vR, p`` (in any order).
Bs : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the bottom leg of the environment, labels are ``vL, vR, p`` (in any order).
R : :class:`~tenpy.linalg.np_conserved.Array`
Initial right environment, labels are ``vL, vR`` (in any order), if 'Ws' is `None`,
otherwise also need a leg ``wL``
Ws : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the middle leg of the environment, labels are ``wL, wR, p, p*`` (in any order).
Returns
-------
temp : :class:`~tenpy.linalg.np_conserved.Array`
The new environment with the tensors above included.
"""
temp = R.copy()
for i in reversed(range(len(As))):
temp = npc.tensordot(Bs[i].conj(), temp, axes=(['vR*'], ['vL*']))
if Ws is not None:
temp = npc.tensordot(Ws[i], temp, axes=(['wR', 'p'], ['wL', 'p*']))
temp = npc.tensordot(As[i], temp, axes=(['vR', 'p'], ['vL', 'p*']))
return temp
def append_left_env(As, Bs, L, Ws=None):
"""Contract all tensors in As and Bs to the left environment L.
If Ws given, contract as a MPO environment.
Parameters
----------
As : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the top leg of the environment, labels are ``vL, vR, p`` (in any order).
Bs : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the bottom leg of the environment, labels are ``vL, vR, p`` (in any order).
L : :class:`~tenpy.linalg.np_conserved.Array`
Initial left environment, labels are ``vL, vR`` (in any order), if 'Ws' is `None`,
otherwise also need a leg ``wR``
Ws : list of :class:`~tenpy.linalg.np_conserved.Array`
Arrays to add to the middle leg of the environment, labels are ``wL, wR, p, p*`` (in any order).
Returns
-------
temp : :class:`~tenpy.linalg.np_conserved.Array`
The new environment with the tensors above included.
"""
temp = L.copy()
for i in range(len(As)):
temp = npc.tensordot(temp, Bs[i].conj(), axes=(['vR*'], ['vL*']))
if Ws is not None:
temp = npc.tensordot(temp, Ws[i], axes=(['wR', 'p*'], ['wL', 'p']))
temp = npc.tensordot(temp, As[i], axes=(['vR', 'p*'], ['vL', 'p']))
return temp
def construct_orthogonal(M, left=True):
"""find (left) orthogonal complement of tensor M
It finds Q such that::
.--Q --
| |
| | = 0
| |
.--M --
if left==False, find the right complement accordingly.
Parameters
----------
M : list of :class:`~tenpy.linalg.np_conserved.Array`
Array for which we want to find the orthogonal complement, labels are ``vL, vR, p`` (in any order).
left : bool
Whether we want to compute the left or right complement.
Returns
-------
Q : :class:`~tenpy.linalg.np_conserved.Array`
The orthogonal complement, such that the relation above is fulfilled.
"""
if left:
M = M.copy().combine_legs([['vL', 'p'], ['vR']], qconj=[+1, -1])
Q = npc.orthogonal_columns(M, 'vR')
assert npc.norm(npc.tensordot(Q, M.conj(), axes=(['(vL.p)'], ['(vL*.p*)']))) < 1.e-12
else:
M = M.copy().combine_legs([['vL'], ['p', 'vR']], qconj=[+1, -1])
Q = npc.orthogonal_columns(M.transpose(['(p.vR)', '(vL)']),
'vL').itranspose(['vL', '(p.vR)'])
assert npc.norm(npc.tensordot(Q, M.conj(), axes=(['(p.vR)'], ['(p*.vR*)']))) < 1.e-12
return Q.split_legs()
class PlaneWaveExcitationEngine(Algorithm):
r""" Base engine to compute quasiparticle excitations for uniform MPS.
Parameters are the same as for :class:`~tenpy.algorithms.algorithm.Algorithm`.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
:include: Algorithm
lanczos_params : dict
Lanczos parameters as described in :cfg:config:`KrylovBased`.
lambda_C1 : float
Energy shift from contracting the infinite environments. If `None`, compute it again.
init_env_data : dict
Dictionary as returned by ``self.env.get_initialization_data()`` from
:meth:`~tenpy.networks.mpo.MPOEnvironment.get_initialization_data`.
Attributes
----------
psi : :class:`~tenpy.networks.uniform_mps.UniformMPS`
The uniform MPS for which we compute (orthogonal) excitations.
model : :class:`~tenpy.models.model.MPOModel`
The model defining the Hamiltonian.
ACs : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'center-site' tensors of psi.
ALs : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'left-orthonormal' tensors of psi.
AR : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'right-orthonormal' tensors of psi.
Cs : list of :class:`~tenpy.linalg.np_conserved.Array`
The center matrices of psi.
Ws : list of :class:`~tenpy.linalg.np_conserved.Array`
The matrices of the MPO.
VLs : list of :class:`~tenpy.linalg.np_conserved.Array`
The orthogonal complements for each AL.
"""
def __init__(self, psi, model, options, **kwargs):
super().__init__(psi, model, options, **kwargs)
assert self.psi.L == self.model.H_MPO.L
self.L = self.psi.L
self.ALs = [self.psi.get_AL(i) for i in range(self.L)]
self.ARs = [self.psi.get_AR(i) for i in range(self.L)]
self.ACs = [self.psi.get_AC(i) for i in range(self.L)]
self.Cs = [self.psi.get_C(i) for i in range(self.L)] # C on the left
self.H = self.model.H_MPO
self.Ws = [self.H.get_W(i) for i in range(self.L)]
if len(self.Ws) < len(self.ALs):
assert len(self.ALs) % len(self.Ws)
self.Ws = self.Ws * len(self.ALs) // len(self.Ws)
self.IdL = self.H.get_IdL(0)
self.IdR = self.H.get_IdR(-1)
self.guess_init_env_data = self.options.get('init_env_data', None)
# Construct VL, needed to parametrize - B - as - VL - X -
# | |
# Use prescription under Eq. 85 in Tangent Space lecture notes.
self.VLs = [construct_orthogonal(self.ALs[i]) for i in range(self.L)]
# Get left and right generalized eigenvalues
self.boundary_env_data, self.energy_density, _ = MPOTransferMatrix.find_init_LP_RP(
self.H, self.psi, calc_E=True, guess_init_env_data=self.guess_init_env_data)
self.energy_density = np.mean(self.energy_density)
self.LW = self.boundary_env_data['init_LP']
self.RW = self.boundary_env_data['init_RP']
# We create GS_env_L and GS_env_R to make topological easier.
self.GS_env = self.GS_env_L = self.GS_env_R = MPOEnvironment(self.psi, self.H, self.psi,
**self.boundary_env_data)
self.lambda_C1 = options.get('lambda_C1', None)
if self.lambda_C1 is None:
C0_L = self.Cs[0]
norm = npc.tensordot(C0_L, C0_L.conj(), axes=(['vL', 'vR'], ['vL*', 'vR*']))
self.lambda_C1 = npc.tensordot(C0_L, self.RW, axes=(['vR'], ['vL']))
self.lambda_C1 = npc.tensordot(self.LW,
self.lambda_C1,
axes=(['wR', 'vR'], ['wL', 'vL']))
self.lambda_C1 = npc.tensordot(
self.lambda_C1, C0_L.conj(), axes=(['vR*', 'vL*'], ['vL*', 'vR*'])) / norm
self.aligned_H = self.Aligned_Effective_H(self)
strange = []
for i in range(self.L):
temp_L = self.GS_env.get_LP(i)
temp_R = self.GS_env.get_RP(i)
temp = append_left_env([self.VLs[i]], [self.ACs[i]], temp_L, Ws=[self.Ws[i]])
temp = npc.tensordot(temp, temp_R, axes=(['wR', 'vR*'], ['wL', 'vL*']))
strange.append(npc.norm(temp))
logger.info("Norm of H|psi> projected into the tangent space on each site: %r.", strange)
def run(self, p, qtotal_change=None, orthogonal_to=[], E_boosts=[], num_ev=1):
""" Run the plane-wave algorithm to find excited states of the given model.
Parameters
----------
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
qtotal_change : list of int
Charge sectors for each of the defined charges.
orthogonal_to : list of list of :class:`~tenpy.linalg.np_conserved.Array`
Find excitations orthogonal to previously found tensors X.
E_boosts: list of float
energy boosts for orthogonal states
num_ev: int
Number of eigenvalues and eigenvectors, that we extract from a single Arnoldi/ Lanczos run
Returns
-------
Es : list of float
Energies of the lowest-energy excitations. Number equal to `num_ev`.
psis : list of :class:`~tenpy.networks.momentum_mps.MomentumMPS`
MomentumMPS corresponding to the lowest-energy excitations. Number equal to `num_ev`.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
E_boost : float
uniform strength of the energy boosts (instead of specifying a list), default to E_boost=100
"""
self.unaligned_H = self.Unaligned_Effective_H(self, p)
effective_H = SumNpcLinearOperator(self.aligned_H, self.unaligned_H)
lanczos_params = self.options.subconfig('lanczos_params')
X_init = self.initial_guess(qtotal_change)
if len(E_boosts) != len(orthogonal_to):
E_boost = self.options.get('E_boost', 100)
E_boosts = [E_boost] * len(orthogonal_to)
if len(orthogonal_to) > 0:
effective_H = BoostNpcLinearOperator(effective_H, E_boosts, orthogonal_to)
if num_ev > 1:
lanczos_params['which'] = 'SR'
lanczos_params['num_ev'] = num_ev
energies, Xs, N = Arnoldi(effective_H, X_init, lanczos_params).run()
psis = []
Es = []
for E, X in zip(energies, Xs):
psis.append(MomentumMPS(X, self.psi, p))
Es.append(E - self.lambda_C1 - self.energy_density * self.L)
else:
energy, X, N = LanczosGroundState(effective_H, X_init, lanczos_params).run()
Es = [energy - self.lambda_C1 - self.energy_density * self.L]
psis = [MomentumMPS(X, self.psi, p)]
if N == lanczos_params.get('N_max', 20):
import warnings
warnings.warn('Maximum Lanczos iterations needed; be wary of results.')
return np.real_if_close(Es), psis, N
def resume_run(self):
raise NotImplementedError()
def energy(self, p, X):
"""
Compute the energy of excited states
Parameters
----------
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
X : list of :class:`~tenpy.linalg.np_conserved.Array`
Excitation tensors for each site of the unit cell.
Returns
-------
energy : float
Energy
"""
self.unaligned_H = self.Unaligned_Effective_H(self, p)
effective_H = SumNpcLinearOperator(self.aligned_H, self.unaligned_H)
HX = effective_H.matvec(X)
E = np.real(npc.inner(X, HX)).item()
return E - self.energy_density * self.L - self.lambda_C1
def infinite_sum_right(self, p, X):
"""
Infinite sum to the right, see Eq. (194) in :cite:`vanderstraeten2019`
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
X : list of :class:`~tenpy.linalg.np_conserved.Array`
Current excitation tensors for each site of the unit cell.
Returns
-------
R_sum : :class:`~tenpy.linalg.np_conserved.Array`
Array representing the right environment including the `B` tensors at each site.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
sum_method : ``explicit`` | ``GMRES``
Whether to explicitly sum the environment by applying the unit cell tensors until convergence (default) or solving the geometric series with the GMRES method.
sum_tol : float
Convergence criterion for the explicit summation.
sum_iterations : int
Maximum number of iterations for the explicit summation (default sum_iterations=100).
"""
sum_tol = self.options.get('sum_tol', 1.e-10)
sum_iterations = self.options.get('sum_iterations', 100)
sum_method = self.options.get('sum_method', 'explicit')
B = npc.tensordot(self.VLs[self.L - 1], X[self.L - 1], axes=(['vR'], ['vL']))
RB = append_right_env([B], [self.ARs[self.L - 1]], self.RW, Ws=[self.Ws[self.L - 1]])
for i in reversed(range(0, self.L - 1)):
B = npc.tensordot(self.VLs[i], X[i], axes=(['vR'], ['vL']))
RB = append_right_env([B], [self.ARs[i]], self.GS_env_R.get_RP(i), Ws=[self.Ws[i]]) + \
append_right_env([self.ALs[i]], [self.ARs[i]], RB, Ws=[self.Ws[i]])
R = RB
if np.isclose(npc.norm(R), 0):
return R
if sum_method == 'explicit':
R_sum = R.copy()
for _ in range(sum_iterations):
R = np.exp(-1.0j * p * self.L) * append_right_env(
self.ALs, self.ARs, R, Ws=self.Ws)
R_sum.iadd_prefactor_other(1., R)
if npc.norm(R) < sum_tol:
break
return R_sum
elif 'GMRES' in sum_method:
class helper_matvec(NpcLinearOperator):
def __init__(self, excit, ALs, ARs, Ws, sum_method):
self.ALs = ALs
self.ARs = ARs
self.Ws = Ws
self.sum_method = sum_method
self.excit = excit
def matvec(self, vec):
Tr = append_right_env(self.ALs, self.ARs, vec, Ws=self.Ws)
if 'reg' in self.sum_method:
raise NotImplementedError(
'GMRES-reg not implemented for multi-site unit cell.')
lr = npc.tensordot(self.excit.l_LR,
vec,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
llr = npc.tensordot(self.excit.LWCc,
vec,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
T1_r = self.excit.r_LR * (
(self.excit.e_LR - 1) * lr + llr) + self.excit.CRW * lr
Tr = Tr - T1_r
return vec - np.exp(-1.0j * p * self.excit.L) * Tr
tm_op = helper_matvec(self, self.ALs, self.ARs, self.Ws, sum_method)
GMRES_params = self.options.subconfig('GMRES_params')
R_sum, _, _, _ = GMRES(tm_op, npc.Array.zeros_like(R) * 1.j, R, GMRES_params).run()
return R_sum
else:
raise ValueError('Sum method', sum_method, 'not recognized!')
def infinite_sum_left(self, p, X):
"""
Infinite sum to the left, see Eq. (194) in :cite:`vanderstraeten2019`
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
X : list of :class:`~tenpy.linalg.np_conserved.Array`
Current excitation tensors for each site of the unit cell.
Returns
-------
L_sum : :class:`~tenpy.linalg.np_conserved.Array`
Array representing the left environment including the `B` tensors at each site.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
sum_method : ``explicit`` | ``GMRES``
Whether to explicitly sum the environment by applying the unit cell tensors until convergence (default) or solving the geometric series with the GMRES method.
sum_tol : float
Convergence criterion for the explicit summation.
sum_iterations : int
Maximum number of iterations for the explicit summation (default sum_iterations=100).
"""
sum_tol = self.options.get('sum_tol', 1.e-10)
sum_iterations = self.options.get('sum_iterations', 100)
sum_method = self.options.get('sum_method', 'explicit')
B = npc.tensordot(self.VLs[0], X[0], axes=(['vR'], ['vL']))
LB = append_left_env([B], [self.ALs[0]], self.LW, Ws=[self.Ws[0]])
for i in range(1, self.L):
B = npc.tensordot(self.VLs[i], X[i], axes=(['vR'], ['vL']))
LB = append_left_env([B], [self.ALs[i]], self.GS_env_L.get_LP(i), Ws=[self.Ws[i]]) + \
append_left_env([self.ARs[i]], [self.ALs[i]], LB, Ws=[self.Ws[i]])
L = LB
if np.isclose(npc.norm(L), 0):
return L
if sum_method == 'explicit':
L_sum = L.copy()
for i in range(sum_iterations):
L = np.exp(1.0j * p * self.L) * append_left_env(self.ARs, self.ALs, L, Ws=self.Ws)
L_sum.iadd_prefactor_other(1., L)
if npc.norm(L) < sum_tol:
break
return L_sum
elif 'GMRES' in sum_method:
class helper_matvec(NpcLinearOperator):
def __init__(self, excit, ALs, ARs, Ws, sum_method):
self.ALs = ALs
self.ARs = ARs
self.Ws = Ws
self.sum_method = sum_method
self.excit = excit
def matvec(self, vec):
lT = append_left_env(self.ARs, self.ALs, vec, Ws=self.Ws)
if 'reg' in self.sum_method:
raise NotImplementedError(
'GMRES-reg not implemented for multi-site unit cell.')
lr = npc.tensordot(vec,
self.excit.r_RL,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
lrr = npc.tensordot(vec,
self.excit.CcRW,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
T1_l = self.excit.l_RL * (
(self.excit.e_RL - 1) * lr + lrr) + self.excit.LWC * lr
lT = lT - T1_l
return vec - np.exp(1.0j * p * self.excit.L) * lT
tm_op = helper_matvec(self, self.ALs, self.ARs, self.Ws, sum_method)
GMRES_params = self.options.subconfig('GMRES_params')
L_sum, _, _, _ = GMRES(tm_op, npc.Array.zeros_like(L) * 1.j, L, GMRES_params).run()
return L_sum
else:
raise ValueError('Sum method', sum_method, 'not recognized!')
class Aligned_Effective_H(NpcLinearOperator):
r"""Class defining the effective Hamiltonian for the excitation tensors `X`.
Where the `B` tensors are in the same unit cell as the tensors we want to update.
For a single-site unit cell the effective Hamiltonian looks like this::
| .--- B ---.
| | | |
| LW----W0----RW
| | | |
| .---VL- --.
Parameters
----------
outer : :class:`PlaneWaveExcitationEngine`
Parent engine for the plane wave excitation ansatz.
"""
def __init__(self, outer):
self.ALs = outer.ALs
self.ARs = outer.ARs
self.VLs = outer.VLs
self.LW = outer.LW
self.RW = outer.RW
self.Ws = outer.Ws
self.outer = outer
def matvec(self, vec):
total_vec = [npc.Array.zeros_like(v) for v in vec]
for i in range(self.outer.L):
LB = npc.Array.zeros_like(self.LW)
RB = npc.Array.zeros_like(self.RW)
for j in range(i):
B = npc.tensordot(self.VLs[j], vec[j], axes=(['vR'], ['vL']))
if j > 0:
LB = append_left_env([B], [self.ALs[j]], self.outer.GS_env_L.get_LP(j), Ws=[self.Ws[j]]) + \
append_left_env([self.ARs[j]], [self.ALs[j]], LB, Ws=[self.Ws[j]]) # Does one extra multiplication when i = 0
else:
LB = append_left_env([B], [self.ALs[j]],
self.outer.GS_env_L.get_LP(j),
Ws=[self.Ws[j]])
B = npc.tensordot(self.VLs[i], vec[i], axes=(['vR'], ['vL']))
LB = append_left_env([self.ARs[i]], [self.VLs[i]], LB, Ws=[self.Ws[i]])
LP1 = append_left_env([self.ALs[i]], [self.VLs[i]],
self.outer.GS_env_L.get_LP(i),
Ws=[self.Ws[i]])
LP2 = append_left_env([B], [self.VLs[i]],
self.outer.GS_env_L.get_LP(i),
Ws=[self.Ws[i]])
for j in reversed(range(i + 1, self.outer.L)):
B = npc.tensordot(self.VLs[j], vec[j], axes=(['vR'], ['vL']))
if j < self.outer.L - 1:
RB = append_right_env([B], [self.ARs[j]], self.outer.GS_env_R.get_RP(j), Ws=[self.Ws[j]]) + \
append_right_env([self.ALs[j]], [self.ARs[j]], RB, Ws=[self.Ws[j]])
else:
RB = append_right_env([B], [self.ARs[j]],
self.outer.GS_env_R.get_RP(j),
Ws=[self.Ws[j]])
if i > 0:
total_vec[i] += npc.tensordot(LB,
self.outer.GS_env_R.get_RP(i),
axes=(['vR', 'wR'], ['vL', 'wL']))
if i < self.outer.L - 1:
total_vec[i] += npc.tensordot(LP1, RB, axes=(['vR', 'wR'], ['vL', 'wL']))
total_vec[i] += npc.tensordot(LP2,
self.outer.GS_env_R.get_RP(i),
axes=(['vR', 'wR'], ['vL', 'wL']))
return total_vec
class Unaligned_Effective_H(NpcLinearOperator):
r"""Class defining the effective Hamiltonian for the excitation tensors `X`, where the `B` tensors are in left (LB) or right (RB) environment.
For a single-site unit cell the effective Hamiltonian looks like this::
| .--- AR ---. .--- AL ---.
| ip | | | -ip | | |
| e LB----W0----RW + e LW----W0----RB
| | | | | | |
| .---VL- --. .---VL- --.
Parameters
----------
outer : :class:`PlaneWaveExcitationEngine`
Parent engine for the plane wave excitation ansatz.
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
"""
def __init__(self, outer, p):
self.ALs = outer.ALs
self.ARs = outer.ARs
self.VLs = outer.VLs
self.LW = outer.LW
self.RW = outer.RW
self.Ws = outer.Ws
self.p = p
self.outer = outer
def matvec(self, vec):
total = [npc.Array.zeros_like(v) for v in vec]
inf_sum_TR = self.outer.infinite_sum_right(self.p, vec)
cached_TR = [inf_sum_TR]
for i in reversed(range(1, self.outer.L)):
cached_TR.insert(
0, append_right_env([self.ALs[i]], [self.ARs[i]],
cached_TR[0],
Ws=[self.Ws[i]]))
for i in range(self.outer.L):
LP_VL = append_left_env([self.ALs[i]], [self.VLs[i]],
self.outer.GS_env_L.get_LP(i),
Ws=[self.Ws[i]])
X_out_left = np.exp(-1.0j * self.p * self.outer.L) * npc.tensordot(
LP_VL, cached_TR[i], axes=(['vR', 'wR'], ['vL', 'wL']))
X_out_left.ireplace_labels(['vR*', 'vL*'], ['vL', 'vR'])
total[i] += X_out_left
cached_TR = []
inf_sum_TL = self.outer.infinite_sum_left(self.p, vec)
cached_TL = [inf_sum_TL]
for i in range(0, self.outer.L - 1):
cached_TL.append(
append_left_env([self.ARs[i]], [self.ALs[i]], cached_TL[-1], Ws=[self.Ws[i]]))
for i in reversed(range(self.outer.L)):
TL_VL = append_left_env([self.ARs[i]], [self.VLs[i]],
cached_TL[i],
Ws=[self.Ws[i]])
X_out_left = np.exp(1.0j * self.p * self.outer.L) * npc.tensordot(
TL_VL, self.outer.GS_env_R.get_RP(i), axes=(['vR', 'wR'], ['vL', 'wL']))
X_out_left.ireplace_labels(['vR*', 'vL*'], ['vL', 'vR'])
total[i] += X_out_left
cached_TL = []
return total
def initial_guess(self, qtotal_change):
"""
Initial guess for the `X` tensors within a fixed charge sector.
Parameters
----------
qtotal_change : list of int
For each charge sector specify how `X` should change the charge.
Returns
-------
X_init : list of :class:`~tenpy.linalg.np_conserved.Array`
Initial guess for excitation tensors for each site of the unit cell.
"""
X_init = []
valid_charge = False
for i in range(self.L):
vL = self.VLs[i].get_leg('vR').conj()
vR = self.ALs[(i + 1) % self.L].get_leg('vL').conj()
th0 = npc.Array.from_func(np.ones, [vL, vR],
dtype=self.psi.dtype,
qtotal=qtotal_change,
labels=['vL', 'vR'])
if np.isclose(npc.norm(th0), 0):
logger.warn("Initial guess for an X is zero; charges not be allowed on site %d.",
i)
else:
valid_charge = True
LP = self.GS_env_L.get_LP(i, store=True)
RP = self.GS_env_R.get_RP(i, store=True)
LP = append_left_env([self.VLs[i]], [self.VLs[i]], LP, Ws=[self.Ws[i]])
H0 = ZeroSiteH.from_LP_RP(LP, RP)
if self.model.H_MPO.explicit_plus_hc:
H0 = SumNpcLinearOperator(H0, H0.adjoint())
lanczos_params = self.options.subconfig('lanczos_params')
_, th0, _ = LanczosGroundState(H0, th0, lanczos_params).run()
X_init.append(th0)
logger.info("Norms of the initial guess: %r.", [npc.norm(x) for x in X_init])
assert valid_charge, "No X is non-zero; charge is not valid for gluing."
return X_init
class MultiSitePlaneWaveExcitationEngine(Algorithm):
r""" Engine to compute quasiparticle excitations across multiple sites for uniform MPS. For each site in the unit cell one multi-site excitation tensor is computed.
Parameters are the same as for :class:`~tenpy.algorithms.algorithm.Algorithm`.
Options
-------
.. cfg:config :: MultiSitePlaneWaveExcitationEngine
:include: Algorithm
lanczos_params : dict
Lanczos parameters as described in :cfg:config:`KrylovBased`.
lambda_C1 : float
Energy shift from contracting the infinite environments. If `None`, compute it again.
init_env_data : dict
Dictionary as returned by ``self.env.get_initialization_data()`` from
:meth:`~tenpy.networks.mpo.MPOEnvironment.get_initialization_data`.
excitation_size : int
Number of sites of the excitation, i.e. how many sites in the uniform MPS are replaced with orthogonal tensors. This can be larger than the unit cell or incommensurate.
Attributes
----------
psi : :class:`~tenpy.networks.uniform_mps.UniformMPS`
The uniform MPS for which we compute (orthogonal) excitations.
model : :class:`~tenpy.models.model.MPOModel`
The model defining the Hamiltonian.
ACs : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'center-site' tensors of psi.
ALs : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'left-orthonormal' tensors of psi.
AR : list of :class:`~tenpy.linalg.np_conserved.Array`
The 'right-orthonormal' tensors of psi.
Cs : list of :class:`~tenpy.linalg.np_conserved.Array`
The center matrices of psi.
Ws : list of :class:`~tenpy.linalg.np_conserved.Array`
The matrices of the MPO.
VLs : list of :class:`~tenpy.linalg.np_conserved.Array`
The orthogonal complements for each AL.
"""
def __init__(self, psi, model, options, **kwargs):
super().__init__(psi, model, options, **kwargs)
assert self.psi.L == self.model.H_MPO.L
self.L = self.psi.L
self.size = self.options.get('excitation_size', 1)
assert self.size >= 1
self.ALs = [self.psi.get_AL(i) for i in range(self.L)]
self.ARs = [self.psi.get_AR(i) for i in range(self.L)]
self.ACs = [self.psi.get_AC(i) for i in range(self.L)]
self.Cs = [self.psi.get_C(i) for i in range(self.L)] # C on the left
self.H = self.model.H_MPO
self.Ws = [self.H.get_W(i) for i in range(self.L)]
if len(self.Ws) < len(self.ALs):
assert len(self.ALs) % len(self.Ws)
self.Ws = self.Ws * len(self.ALs) // len(self.Ws)
self.IdL = self.H.get_IdL(0)
self.IdR = self.H.get_IdR(-1)
self.guess_init_env_data = self.options.get('init_env_data', None)
# Construct VL, needed to parametrize - B - as - VL - X -
# | |
# Use prescription under Eq. 85 in Tangent Space lecture notes.
self.VLs = [construct_orthogonal(self.ALs[i]) for i in range(self.L)]
# Get left and right generalized eigenvalues
self.boundary_env_data, self.energy_density, _ = MPOTransferMatrix.find_init_LP_RP(
self.H, self.psi, calc_E=True, guess_init_env_data=self.guess_init_env_data)
self.energy_density = np.mean(self.energy_density)
self.LW = self.boundary_env_data['init_LP']
self.RW = self.boundary_env_data['init_RP']
# We create GS_env_L and GS_env_R to make topological easier.
self.GS_env = self.GS_env_L = self.GS_env_R = MPOEnvironment(self.psi, self.H, self.psi,
**self.boundary_env_data)
self.lambda_C1 = options.get('lambda_C1', None)
if self.lambda_C1 is None:
C0_L = self.Cs[0]
norm = npc.tensordot(C0_L, C0_L.conj(), axes=(['vL', 'vR'], ['vL*', 'vR*']))
self.lambda_C1 = npc.tensordot(C0_L, self.RW, axes=(['vR'], ['vL']))
self.lambda_C1 = npc.tensordot(self.LW,
self.lambda_C1,
axes=(['wR', 'vR'], ['wL', 'vL']))
self.lambda_C1 = npc.tensordot(
self.lambda_C1, C0_L.conj(), axes=(['vR*', 'vL*'], ['vL*', 'vR*'])) / norm
strange = []
for i in range(self.L):
temp_L = self.GS_env.get_LP(i)
temp_R = self.GS_env.get_RP(i)
temp = append_left_env([self.VLs[i]], [self.ACs[i]], temp_L, Ws=[self.Ws[i]])
temp = npc.tensordot(temp, temp_R, axes=(['wR', 'vR*'], ['wL', 'vL*']))
strange.append(npc.norm(temp))
logger.info("Norm of H|psi> projected into the tangent space on each site: %r.", strange)
def run(self, p, qtotal_change=None, orthogonal_to=[], E_boosts=[], num_ev=1):
""" Run the plane-wave algorithm to find excited states of the given model.
Parameters
----------
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
qtotal_change : list of int
Charge sectors for each of the defined charges.
orthogonal_to : list of list of :class:`~tenpy.linalg.np_conserved.Array`
Find excitations orthogonal to previously found tensors X.
E_boosts: list of float
energy boosts for orthogonal states
num_ev: int
Number of eigenvalues and eigenvectors, that we extract from a single Arnoldi/ Lanczos run
Returns
-------
Es : list of float
Energies of the lowest-energy excitations. Number equal to `num_ev`.
psis : list of :class:`~tenpy.networks.momentum_mps.MomentumMPS`
MomentumMPS corresponding to the lowest-energy excitations. Number equal to `num_ev`.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
E_boost : float
uniform strength of the energy boosts (instead of specifying a list), default to E_boost=100
"""
self.aligned_H = self.Aligned_Effective_H(self, p)
self.unaligned_H = self.Unaligned_Effective_H(self, p)
effective_H = SumNpcLinearOperator(self.aligned_H, self.unaligned_H)
lanczos_params = self.options.subconfig('lanczos_params')
X_init = self.initial_guess(qtotal_change)
if len(E_boosts) != len(orthogonal_to):
E_boost = self.options.get('E_boost', 100)
E_boosts = [E_boost] * len(orthogonal_to)
if len(orthogonal_to) > 0:
effective_H = BoostNpcLinearOperator(effective_H, E_boosts, orthogonal_to)
# how many unit cells to include all excitations
multiple_unit_cell = int(np.ceil((self.L - 1 + self.size) / self.L))
if num_ev > 1:
lanczos_params['which'] = 'SR'
lanczos_params['num_ev'] = num_ev
energies, Xs, N = Arnoldi(effective_H, X_init, lanczos_params).run()
psis = []
Es = []
for E, X in zip(energies, Xs):
psis.append(MomentumMPS(X, self.psi, p, self.size))
Es.append(E - self.lambda_C1 - self.energy_density * (self.L * multiple_unit_cell))
else:
energy, X, N = LanczosGroundState(effective_H, X_init, lanczos_params).run()
Es = [energy - self.lambda_C1 - self.energy_density * (self.L * multiple_unit_cell)]
psis = [MomentumMPS(X, self.psi, p, self.size)]
if N == lanczos_params.get('N_max', 20):
import warnings
warnings.warn('Maximum Lanczos iterations needed; be wary of results.')
return np.real_if_close(Es), psis, N
def resume_run(self):
raise NotImplementedError()
def energy(self, p, X):
"""
Compute the energy of excited states
Parameters
----------
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
X : list of :class:`~tenpy.linalg.np_conserved.Array`
Excitation tensors for each site of the unit cell.
Returns
-------
energy : float
Energy
"""
multiple_unit_cell = int(np.ceil((self.L - 1 + self.size) / self.L))
self.aligned_H = self.Aligned_Effective_H(self, p)
self.unaligned_H = self.Unaligned_Effective_H(self, p)
effective_H = SumNpcLinearOperator(self.aligned_H, self.unaligned_H)
HX = effective_H.matvec(X)
E = np.real(npc.inner(X, HX)).item()
return E - self.lambda_C1 - self.energy_density * (self.L * multiple_unit_cell)
def attach_right(self, VL, X, As, R, Ws=None):
"""
attach excitation tensors to a right environment
"""
B = npc.tensordot(VL.replace_label('p', 'p0'), X, axes=(['vR'], ['vL']))
RB = npc.tensordot(B, R, axes=(['vR'], ['vL']))
for i in reversed(range(len(As))):
p = 'p' + str(i)
if Ws is not None:
RB = npc.tensordot(RB, Ws[i], axes=([p, 'wL'], ['p*', 'wR']))
RB = npc.tensordot(RB, As[i].conj(), axes=(['p', 'vL*'], ['p*', 'vR*']))
return RB
def _starting_right_TR(self, X):
"""
Sum up all single-B environments from the right and fill with tensors to complete unit cell
"""
i = 0
RP = self.GS_env_R.get_RP(i + self.size - 1)
RB = self.attach_right(self.VLs[i],
X[i], [self.ARs[j % self.L] for j in range(i, i + self.size)],
RP,
Ws=[self.Ws[j % self.L] for j in range(i, i + self.size)])
RB = append_right_env(self.ALs[:i], self.ARs[:i], RB, Ws=self.Ws[:i])
RW = RB
for i in range(1, self.L):
RP = self.GS_env_R.get_RP(i + self.size - 1)
RB = self.attach_right(self.VLs[i],
X[i], [self.ARs[j % self.L] for j in range(i, i + self.size)],
RP,
Ws=[self.Ws[j % self.L] for j in range(i, i + self.size)])
RB = append_right_env(self.ALs[:i], self.ARs[:i], RB, Ws=self.Ws[:i])
RW += RB
return RW
def infinite_sum_right(self, p, X):
"""
Infinite sum to the right, see Eq. (194) in :cite:`vanderstraeten2019`
p : float
The momentum of the state; for unit cells larger than 1, we already include the
factor of the smaller Brillouin zone: p*L.
X : list of :class:`~tenpy.linalg.np_conserved.Array`
Current excitation tensors for each site of the unit cell.
Returns
-------
R_sum : :class:`~tenpy.linalg.np_conserved.Array`
Array representing the right environment including the `B` tensors at each site.
Options
-------
.. cfg:config :: PlaneWaveExcitationEngine
sum_method : ``explicit`` | ``GMRES``
Whether to explicitly sum the environment by applying the unit cell tensors until convergence (default) or solving the geometric series with the GMRES method.
sum_tol : float
Convergence criterion for the explicit summation.
sum_iterations : int
Maximum number of iterations for the explicit summation (default sum_iterations=100).
"""
sum_tol = self.options.get('sum_tol', 1.e-10)
sum_iterations = self.options.get('sum_iterations', 100)
sum_method = self.options.get('sum_method', 'explicit')
R = self._starting_right_TR(X)
if np.isclose(npc.norm(R), 0):
return R
if sum_method == 'explicit':
R_sum = R.copy()
for _ in range(sum_iterations):
R = np.exp(-1.0j * p * self.L) * append_right_env(
self.ALs, self.ARs, R, Ws=self.Ws)
R_sum.iadd_prefactor_other(1., R)
if npc.norm(R) < sum_tol:
break
return R_sum
elif 'GMRES' in sum_method:
class helper_matvec(NpcLinearOperator):
def __init__(self, excit, ALs, ARs, Ws, sum_method):
self.ALs = ALs
self.ARs = ARs
self.Ws = Ws
self.sum_method = sum_method
self.excit = excit
def matvec(self, vec):
Tr = append_right_env(self.ALs, self.ARs, vec, Ws=self.Ws)
if 'reg' in self.sum_method:
raise NotImplementedError(
'GMRES-reg not implemented for multi-site unit cell.')
lr = npc.tensordot(self.excit.l_LR,
vec,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
llr = npc.tensordot(self.excit.LWCc,
vec,
axes=(['vR', 'wR', 'vR*'], ['vL', 'wL', 'vL*']))
T1_r = self.excit.r_LR * (
(self.excit.e_LR - 1) * lr + llr) + self.excit.CRW * lr
Tr = Tr - T1_r
return vec - np.exp(-1.0j * p * self.excit.L) * Tr
tm_op = helper_matvec(self, self.ALs, self.ARs, self.Ws, sum_method)
GMRES_params = self.options.subconfig('GMRES_params')
R_sum, _, _, _ = GMRES(tm_op, npc.Array.zeros_like(R) * 1.j, R, GMRES_params).run()
return R_sum
else:
raise ValueError('Sum method', sum_method, 'not recognized!')
def attach_left(self, VL, X, As, L, Ws=None):
"""
attach excitation tensors to a left environment
"""
B = npc.tensordot(VL.replace_label('p', 'p0'), X, axes=(['vR'], ['vL']))
LB = npc.tensordot(L, B, axes=(['vR'], ['vL']))
for i in range(len(As)):