/
dmrg.py
1931 lines (1726 loc) · 86.1 KB
/
dmrg.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
"""Density Matrix Renormalization Group (DMRG).
Although it was originally not formulated with tensor networks,
the DMRG algorithm (invented by Steven White in 1992 :cite:`white1992`) opened the whole field
with its enormous success in finding ground states in 1D.
We implement DMRG in the modern formulation of matrix product states :cite:`schollwoeck2011`,
both for finite systems (``'finite'`` or ``'segment'`` boundary conditions)
and in the thermodynamic limit (``'infinite'`` b.c.).
The function :func:`run` - well - runs one DMRG simulation.
Internally, it generates an instance of an :class:`Sweep`.
This class implements the common functionality like defining a `sweep`,
but leaves the details of the contractions to be performed to the derived classes.
Currently, there are two derived classes implementing the contractions: :class:`SingleSiteDMRGEngine`
and :class:`TwoSiteDMRGEngine`. They differ (as their name implies) in the number of sites which
are optimized simultaneously.
They should both give the same results (up to rounding errors). However, if started from a product
state, :class:`SingleSiteDMRGEngine` depends critically on the use of a :class:`Mixer`, while
:class:`TwoSiteDMRGEngine` is in principle more computationally expensive to run and has
occasionally displayed some convergence issues..
Which one is preffered in the end is not obvious a priori and might depend on the used model.
Just try both of them.
A :class:`Mixer` should be used initially to avoid that the algorithm gets stuck in local energy
minima, and then slowly turned off in the end. For :class:`SingleSiteDMRGEngine`, using a mixer is
crucial, as the one-site algorithm cannot increase the MPS bond dimension by itself.
A generic protocol for approaching a physics question using DMRG is given in :doc:`/intro/protocol`.
"""
# Copyright 2018-2021 TeNPy Developers, GNU GPLv3
import numpy as np
import time
import warnings
import logging
logger = logging.getLogger(__name__)
from ..linalg import np_conserved as npc
from ..networks.mps import MPSEnvironment
from ..networks.mpo import MPOEnvironment
from ..linalg.lanczos import lanczos, lanczos_arpack
from .truncation import truncate, svd_theta
from ..tools.params import asConfig
from ..tools.process import memory_usage
from .mps_common import Sweep, OneSiteH, TwoSiteH
__all__ = [
'run', 'DMRGEngine', 'SingleSiteDMRGEngine', 'TwoSiteDMRGEngine', 'EngineCombine',
'EngineFracture', 'Mixer', 'SingleSiteMixer', 'TwoSiteMixer', 'DensityMatrixMixer', 'chi_list',
'full_diag_effH'
]
def run(psi, model, options):
r"""Run the DMRG algorithm to find the ground state of the given model.
Parameters
----------
psi : :class:`~tenpy.networks.mps.MPS`
Initial guess for the ground state, which is to be optimized in-place.
model : :class:`~tenpy.models.MPOModel`
The model representing the Hamiltonian for which we want to find the ground state.
options : dict
Further optional parameters as described in :cfg:config:`DMRGEngine`.
Returns
-------
info : dict
A dictionary with keys ``'E', 'shelve', 'bond_statistics', 'sweep_statistics'``
Options
-------
.. cfg:config :: DMRG
:include: SingleSiteDMRGEngine, TwoSiteDMRGEngine
active_sites : 1 | 2
The number of active sites to be used by DMRG.
If set to 1, :class:`SingleSiteDMRGEngine` is used.
If set to 2, DMRG is handled by :class:`TwoSiteDMRGEngine`.
"""
# initialize the engine
options = asConfig(options, 'DMRG')
active_sites = options.get('active_sites', 2)
if active_sites == 1:
engine = SingleSiteDMRGEngine(psi, model, options)
elif active_sites == 2:
engine = TwoSiteDMRGEngine(psi, model, options)
else:
raise ValueError("For DMRG, can only use 1 or 2 active sites, not {}".format(active_sites))
E, _ = engine.run()
return {
'E': E,
'shelve': engine.shelve,
'bond_statistics': engine.update_stats,
'sweep_statistics': engine.sweep_stats
}
class Mixer:
"""Base class of a general Mixer.
Since DMRG performs only local updates of the state, it can get stuck in "local minima",
in particular if the Hamiltonain is long-range -- which is the case if one
maps a 2D system ("infinite cylinder") to 1D -- or if one wants to do single-site updates
(currently not implemented in TeNPy).
The idea of the mixer is to perturb the state with the terms of the Hamiltonian
which have contributions in both the "left" and "right" side of the system.
In that way, it adds fluctuation of the quantum numbers and non-zero contributions of the
long-range terms - leading to a significantly improved convergence of DMRG.
The strength of the perturbation is given by the `amplitude` of the mixer.
A good strategy is to choose an initially significant amplitude and let it decay until
the perturbation becomes completely irrelevant and the mixer gets disabled.
This original idea of the mixer was introduced in :cite:`white2005`.
:cite:`hubig2015` discusses the mixer and provides an improved version.
Parameters
----------
env : :class:`~tenpy.networks.mpo.MPOEnvironment`
Environment for contraction ``<psi|H|psi>`` for later
options : dict
Optional parameters as described in the following table.
see :cfg:config:`Mixer`
Options
-------
.. cfg:config :: Mixer
amplitude : float
Initial strength of the mixer. (Should be << 1.)
decay : float
To slowly turn off the mixer, we divide `amplitude` by `decay`
after each sweep. (Should be >= 1.)
disable_after : int
We disable the mixer completely after this number of sweeps.
Attributes
----------
amplitude : float
Current amplitude for mixing.
decay : float
Factor by which `amplitude` is divided after each sweep.
disable_after : int
The number of sweeps after which the mixer should be disabled.
"""
def __init__(self, options):
self.options = options = asConfig(options, 'Mixer')
self.amplitude = options.get('amplitude', 1.e-5)
assert self.amplitude <= 1.
self.decay = options.get('decay', 2.)
assert self.decay >= 1.
if self.decay == 1.:
warnings.warn("Mixer with decay=1. doesn't decay")
self.disable_after = options.get('disable_after', 15)
def update_amplitude(self, sweeps):
"""Update the amplitude, possibly disable the mixer.
Parameters
----------
sweeps : int
The number of performed sweeps, to check if we need to disable the mixer.
Returns
-------
mixer : :class:`Mixer` | None
Returns `self` if we should continue mixing, or ``None``, if the mixer
should be disabled.
"""
self.amplitude /= self.decay
if sweeps >= self.disable_after or self.amplitude <= np.finfo('float').eps:
logger.info("disable mixer after %(sweeps)d sweeps, final amplitude %(amp).2e", {
'sweeps': sweeps,
'amp': self.amplitude
})
return None # disable mixer
return self
def perturb_svd(self, engine, theta, i0, update_LP, update_RP):
"""Perturb the wave function and perform an SVD with truncation.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
The optimized wave function, prepared for svd.
i0 : int
Site index; `theta` lives on ``i0, i0+1``.
update_LP : bool
Whether to calculate the next ``env.LP[i0+1]``.
update_RP : bool
Whether to calculate the next ``env.RP[i0]``.
Returns
-------
U : :class:`~tenpy.linalg.np_conserved.Array`
Left-canonical part of `theta`. Labels ``'(vL.p0)', 'vR'``.
S : 1D ndarray | 2D :class:`~tenpy.linalg.np_conserved.Array`
Without mixer just the singluar values of the array; with mixer it might be a general
matrix; see comment above.
VH : :class:`~tenpy.linalg.np_conserved.Array`
Right-canonical part of `theta`. Labels ``'vL', '(vR.p1)'``.
err : :class:`~tenpy.algorithms.truncation.TruncationError`
The truncation error introduced.
"""
raise NotImplementedError("This function should be implemented in derived classes")
class SingleSiteMixer(Mixer):
"""Mixer for single-site DMRG.
Performs a subspace expansion following :cite:`hubig2015`.
"""
def perturb_svd(self, engine, theta, i0, move_right, next_B):
"""Mix extra terms to theta and perform an SVD.
We calculate the left and right reduced density matrix using the mixer
(which might include applications of `H`).
These density matrices are diagonalized and truncated such that we effectively perform
a svd for the case ``mixer.amplitude=0``.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
The optimized wave function, prepared for svd.
i0 : int
The site index where `theta` lives.
move_right : bool
Whether we move to the right (``True``) or left (``False``).
next_B : :class:`~tenpy.linalg.np_conserved.Array`
The subspace expansion requires to change the tensor on the next site as well.
If `move_right`, it should correspond to ``engine.psi.get_B(i0+1, form='B')``.
If not `move_right`, it should correspond to ``engine.psi.get_B(i0-1, form='A')``.
Returns
-------
U : :class:`~tenpy.linalg.np_conserved.Array`
Left-canonical part of `tensordot(theta, next_B)`. Labels ``'(vL.p0)', 'vR'``.
S : 1D ndarray
(Perturbed) singular values on the new bond (between `theta` and `next_B`).
VH : :class:`~tenpy.linalg.np_conserved.Array`
Right-canonical part of `tensordot(theta, next_B)`. Labels ``'vL', '(p1.vR)'``.
err : :class:`~tenpy.algorithms.truncation.TruncationError`
The truncation error introduced.
"""
theta, next_B = self.subspace_expand(engine, theta, i0, move_right, next_B)
qtotal_LR = [theta.qtotal, None] if move_right else [None, theta.qtotal]
U, S, VH, err, _ = svd_theta(theta,
engine.trunc_params,
qtotal_LR=qtotal_LR,
inner_labels=['vR', 'vL'])
if move_right:
VH = npc.tensordot(VH, next_B, axes=['vR', 'vL'])
else:
U = npc.tensordot(next_B, U, axes=['vR', 'vL'])
return U, S, VH, err
def subspace_expand(self, engine, theta, i0, move_right, next_B):
"""Expand the MPS subspace, to allow the bond dimension to increase.
This is the subspace expansion following :cite:`hubig2015`.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
Optimized guess for the ground state of the effective local Hamiltonian.
i0 : int
Site index at which the local update has taken place.
move_right : bool
Whether the next `i0` of the sweep will be right or left of the current one.
next_B : :class:`~tenpy.linalg.np_conserved.Array`
The subspace expansion requires to change the tensor on the next site as well.
If `move_right`, it should correspond to ``engine.psi.get_B(i0+1, form='B')``.
If not `move_right`, it should correspond to ``engine.psi.get_B(i0-1, form='A')``.
Returns
-------
theta :
Local MPS tensor at site `i0` after subspace expansion.
next_B :
MPS tensor at site `i0+1` or `i0-1` (depending on sweep direction) after subspace
expansion.
"""
eff_H = engine.eff_H
if not engine.combine: # Need to get Heff's even if combine=False
eff_H.combine_Heff()
if move_right: # theta has legs (vL.p0), vR
LHeff = eff_H.LHeff
expand = npc.tensordot(LHeff, theta, axes=['(vR.p0*)', '(vL.p0)'])
expand = expand.combine_legs(['wR', 'vR'], qconj=-1, new_axes=1)
expand *= self.amplitude
theta = npc.concatenate([theta, expand], axis=1, copy=False)
next_B = next_B.extend('vL', expand.legs[1].conj())
else: # theta has legs vL, (p0.vR)
RHeff = eff_H.RHeff
expand = npc.tensordot(theta, RHeff, axes=['(p0.vR)', '(p0*.vL)'])
expand = expand.combine_legs(['vL', 'wL'], qconj=+1)
expand *= self.amplitude
theta = npc.concatenate([theta, expand], axis=0, copy=False)
next_B = next_B.extend('vR', expand.legs[0].conj())
return theta, next_B
class TwoSiteMixer(SingleSiteMixer):
"""Mixer for two-site DMRG.
This is the two-site version of the mixer described in :cite:`hubig2015`.
Equivalent to the :class:`DensityMatrixMixer`, but never construct the full density matrix.
.. todo :
This is still under development.
Works only with :class:`TwoSiteDMRGEngine`.
Has not been ported to `Sweep`-based setup yet. Do we need to?
"""
def perturb_svd(self, engine, theta, i0, move_right):
"""Mix extra terms to theta and perform an SVD.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
The optimized wave function, prepared for svd.
i0 : int
Site index; `theta` lives on ``i0, i0+1``.
update_LP : bool
Whether to calculate the next ``env.LP[i0+1]``.
update_RP : bool
Whether to calculate the next ``env.RP[i0]``.
Returns
-------
U : :class:`~tenpy.linalg.np_conserved.Array`
Left-canonical part of `theta`. Labels ``'(vL.p0)', 'vR'``.
S : 1D ndarray | 2D :class:`~tenpy.linalg.np_conserved.Array`
Without mixer just the singluar values of the array; with mixer it might be a general
matrix; see comment above.
VH : :class:`~tenpy.linalg.np_conserved.Array`
Right-canonical part of `theta`. Labels ``'vL', '(vR.p1)'``.
err : :class:`~tenpy.algorithms.truncation.TruncationError`
The truncation error introduced.
"""
# first perform an SVD as if the mixer didn't exist
qtotal_i0 = engine.psi.get_B(i0, form=None).qtotal
U, S, VH, err, _ = svd_theta(theta,
engine.trunc_params,
qtotal_LR=[qtotal_i0, None],
inner_labels=['vR', 'vL'])
if move_right: # move to the right
U, S, VH, err2 = SingleSiteMixer.perturb_svd(self, engine, U.iscale_axis(S, 1), i0,
move_right, VH)
else: # update_RP is True
U, S, VH, err2 = SingleSiteMixer.perturb_svd(self, engine, VH.iscale_axis(S, 0), i0,
move_right, U)
return U, S, VH, err + err2
class DensityMatrixMixer(Mixer):
"""Mixer based on density matrices.
This mixer constructs density matrices as described in the original paper :cite:`white2005`.
"""
def perturb_svd(self, engine, theta, i0, update_LP, update_RP):
"""Mix extra terms to theta and perform an SVD.
We calculate the left and right reduced density using the mixer
(which might include applications of `H`).
These density matrices are diagonalized and truncated such that we effectively perform
a svd for the case ``mixer.amplitude=0``.
Parameters
----------
engine : :class:`SingleSiteDMRGEngine` | :class:`TwoSiteDMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
The optimized wave function, prepared for svd.
i0 : int
Site index; `theta` lives on ``i0, i0+1``.
update_LP : bool
Whether to calculate the next ``env.LP[i0+1]``.
update_RP : bool
Whether to calculate the next ``env.RP[i0]``.
Returns
-------
U : :class:`~tenpy.linalg.np_conserved.Array`
Left-canonical part of `theta`. Labels ``'(vL.p0)', 'vR'``.
S : 1D ndarray | 2D :class:`~tenpy.linalg.np_conserved.Array`
Without mixer just the singluar values of the array; with mixer it might be a general
matrix; see comment above.
VH : :class:`~tenpy.linalg.np_conserved.Array`
Right-canonical part of `theta`. Labels ``'vL', '(p1.vR)'``.
err : :class:`~tenpy.algorithms.truncation.TruncationError`
The truncation error introduced.
"""
rho_L = self.mix_rho_L(engine, theta, i0, update_LP)
# don't mix left parts, when we're going to the right
rho_L.itranspose(['(vL.p0)', '(vL*.p0*)']) # just to be sure of the order
rho_R = self.mix_rho_R(engine, theta, i0, update_RP)
rho_R.itranspose(['(p1.vR)', '(p1*.vR*)']) # just to be sure of the order
# consider the SVD `theta = U S V^H` (with real, diagonal S>0)
# rho_L ~= theta theta^H = U S V^H V S U^H = U S S U^H (for mixer -> 0)
# Thus, rho_L U = U S S, i.e. columns of U are the eigenvectors of rho_L,
# eigenvalues are S^2.
val_L, U = npc.eigh(rho_L)
U.legs[1] = U.legs[1].to_LegCharge() # explicit conversion: avoid warning in `iproject`
U.iset_leg_labels(['(vL.p0)', 'vR'])
val_L[val_L < 0.] = 0. # for stability reasons
val_L /= np.sum(val_L)
keep_L, _, errL = truncate(np.sqrt(val_L), engine.trunc_params)
U.iproject(keep_L, axes='vR') # in place
U = U.gauge_total_charge(1, engine.psi.get_B(i0, form=None).qtotal)
# rho_R ~= theta^T theta^* = V^* S U^T U* S V^T = V^* S S V^T (for mixer -> 0)
# Thus, rho_L V^* = V^* S S, i.e. columns of V^* are eigenvectors of rho_L
val_R, Vc = npc.eigh(rho_R)
Vc.legs[1] = Vc.legs[1].to_LegCharge()
Vc.iset_leg_labels(['(p1.vR)', 'vL'])
VH = Vc.itranspose(['vL', '(p1.vR)'])
val_R[val_R < 0.] = 0. # for stability reasons
val_R /= np.sum(val_R)
keep_R, _, err_R = truncate(np.sqrt(val_R), engine.trunc_params)
VH.iproject(keep_R, axes='vL')
VH = VH.gauge_total_charge(0, engine.psi.get_B(i0 + 1, form=None).qtotal)
# calculate S = U^H theta V
theta = npc.tensordot(U.conj(), theta, axes=['(vL*.p0*)', '(vL.p0)']) # axes 0, 0
theta = npc.tensordot(theta, VH.conj(), axes=['(p1.vR)', '(p1*.vR*)']) # axes 1, 1
theta.ireplace_labels(['vR*', 'vL*'], ['vL', 'vR']) # for left/right
# normalize `S` (as in svd_theta) to avoid blowing up numbers
theta /= np.linalg.norm(npc.svd(theta, compute_uv=False))
return U, theta, VH, errL + err_R
def mix_rho_L(self, engine, theta, i0, mix_enabled):
"""Calculated mixed reduced density matrix for left site.
Pictorially::
| mix_enabled=False mix_enabled=True
|
| .---theta---. .---theta-------.
| | | | | | | | |
| | | LP---H0--H1--. |
| | | | | | | | | |
| .---theta*--. | xR |
| | | | | |
| LP*--H0*-H1*-. |
| | | | |
| .---theta*------.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
Ground state of the effective Hamiltonian, prepared for svd.
i0 : int
Site index; `theta` lives on ``i0, i0+1``.
mix_enabled : bool
Whether we should perturb the density matrix.
Returns
-------
rho_L : :class:`~tenpy.linalg.np_conserved.Array`
A (hermitian) square array with labels ``'(vL.p0)', '(vL*.p0*)'``,
Mainly the reduced density matrix of the left part, but with some additional mixing.
"""
if not mix_enabled:
return npc.tensordot(theta, theta.conj(), axes=['(p1.vR)', '(p1*.vR*)'])
H = engine.env.H
try:
LHeff = engine.LHeff
except AttributeError:
# TODO: needed?
H0 = H.get_W(i0).replace_labels(['p', 'p*'], ['p0', 'p0*'])
LP = engine.env.get_LP(i0, store=False)
LHeff = npc.tensordot(LP, H0, axes=['wR', 'wL'])
pipeL = theta.get_leg('(vL.p0)')
LHeff = LHeff.combine_legs([['vR*', 'p0'], ['vR', 'p0*']],
pipes=[pipeL, pipeL.conj()],
new_axes=[0, -1])
rho = npc.tensordot(LHeff, theta, axes=['(vR.p0*)', '(vL.p0)']).split_legs('(p1.vR)')
rho_c = rho.conj()
H1 = H.get_W(i0 + 1).replace_labels(['p', 'p*'], ['p1', 'p1*'])
mixer_xR, add_separate_Id = self.get_xR(H1.get_leg('wR'), H.get_IdL(i0 + 2),
H.get_IdR(i0 + 1))
H1m = npc.tensordot(H1, mixer_xR, axes=['wR', 'wL'])
H1m = npc.tensordot(H1m, H1.conj(), axes=[['p1', 'wL*'], ['p1*', 'wR*']])
rho = npc.tensordot(rho, H1m, axes=[['wR', 'p1'], ['wL', 'p1*']])
rho = npc.tensordot(rho, rho_c, axes=(['p1', 'wL*', 'vR'], ['p1*', 'wR*', 'vR*']))
rho.ireplace_labels(['(vR*.p0)', '(vR.p0*)'], ['(vL.p0)', '(vL*.p0*)'])
if add_separate_Id:
rho = rho + npc.tensordot(theta, theta.conj(), axes=['(p1.vR)', '(p1*.vR*)'])
return rho
def mix_rho_R(self, engine, theta, i0, mix_enabled):
"""Calculated mixed reduced density matrix for left site.
Pictorially::
| mix_enabled=False mix_enabled=True
|
| .---theta---. .------theta---.
| | | | | | | | |
| | | | .--H0--H1--RP
| | | | | | | | | |
| .---theta*--. | wL |
| | | | | |
| | .--H0*-H1*-RP*
| | | | |
| .------theta*--.
Parameters
----------
engine : :class:`DMRGEngine`
The DMRG engine calling the mixer.
theta : :class:`~tenpy.linalg.np_conserved.Array`
Ground state of the effective Hamiltonian, prepared for svd.
i0 : int
Site index; `theta` lives on ``i0, i0+1``.
mix_enabled : bool
Whether we should perturb the density matrix.
Returns
-------
rho_R : :class:`~tenpy.linalg.np_conserved.Array`
A (hermitian) square array with labels ``'(p1.vR)', '(p1*.vR*)'``.
Mainly the reduced density matrix of the right part, but with some additional mixing.
"""
if not mix_enabled:
return npc.tensordot(theta, theta.conj(), axes=[['(vL.p0)'], ['(vL*.p0*)']])
H = engine.env.H
try:
RHeff = engine.RHeff
except AttributeError:
H1 = H.get_W(i0 + 1).replace_labels(['p', 'p*'], ['p1', 'p1*'])
RP = engine.env.get_RP(i0 + 1, store=False)
RHeff = npc.tensordot(RP, H1, axes=['wL', 'wR'])
pipeR = theta.get_leg('(p1.vR)')
RHeff = RHeff.combine_legs([['p1', 'vL*'], ['p1*', 'vL']],
pipes=[pipeR, pipeR.conj()],
new_axes=[-1, 0])
rho = npc.tensordot(RHeff, theta, axes=['(p1*.vL)', '(p1.vR)']).split_legs('(vL.p0)')
rho_c = rho.conj()
H0 = H.get_W(i0).replace_labels(['p', 'p*'], ['p0', 'p0*'])
mixer_xL, add_separate_Id = self.get_xL(H0.get_leg('wL'), H.get_IdL(i0), H.get_IdR(i0 - 1))
H0m = npc.tensordot(mixer_xL, H0, axes=['wR', 'wL'])
H0m = npc.tensordot(H0m, H0.conj(), axes=[['wR*', 'p0'], ['wL*', 'p0*']])
rho = npc.tensordot(H0m, rho, axes=[['p0*', 'wR'], ['p0', 'wL']])
rho = npc.tensordot(rho, rho_c, axes=(['p0', 'wR*', 'vL'], ['p0*', 'wL*', 'vL*']))
rho.ireplace_labels(['(p1.vL*)', '(p1*.vL)'], ['(p1.vR)', '(p1*.vR*)'])
if add_separate_Id:
rho = rho + npc.tensordot(theta, theta.conj(), axes=['(vL.p0)', '(vL*.p0*)'])
return rho
def get_xR(self, wR_leg, Id_L, Id_R):
"""Generate the coupling of the MPO legs for the reduced density matrix.
Parameters
----------
wR_leg : :class:`~tenpy.linalg.charges.LegCharge`
LegCharge to be connected to.
IdL : int | ``None``
Index within the leg for which the MPO has only identities to the left.
IdR : int | ``None``
Index within the leg for which the MPO has only identities to the right.
Returns
-------
mixed_xR : :class:`~tenpy.linalg.np_conserved.Array`
Connection of the MPOs on the right for the reduced density matrix `rhoL`.
Labels ``('wL', 'wL*')``.
add_separate_Id : bool
If Id_L is ``None``, we can't include the identity into `mixed_xR`,
so it has to be added directly in :meth:`mix_rho_L`.
"""
x = self.amplitude * np.ones(wR_leg.ind_len, dtype=float)
separate_Id = Id_L is None
if not separate_Id:
x[Id_L] = 1.
if Id_R is not None:
x[Id_R] = 0.
x = npc.diag(x, wR_leg, labels=['wL*', 'wL'])
return x, separate_Id
def get_xL(self, wL_leg, Id_L, Id_R):
"""Generate the coupling of the MPO legs for the reduced density matrix.
Parameters
----------
wL_leg : :class:`~tenpy.linalg.charges.LegCharge`
LegCharge to be connected to.
Id_L : int | ``None``
Index within the leg for which the MPO has only identities to the left.
Id_R : int | ``None``
Index within the leg for which the MPO has only identities to the right.
Returns
-------
mixed_xL : :class:`~tenpy.linalg.np_conserved.Array`
Connection of the MPOs on the left for the reduced density matrix `rhoR`.
Labels ``('wR', 'wR*')``.
add_separate_Id : bool
If Id_R is ``None``, we can't include the identity into `mixed_xL`,
so it has to be added directly in :meth:`mix_rho_R`.
"""
x = self.amplitude * np.ones(wL_leg.ind_len, dtype=float)
separate_Id = Id_R is None
if not separate_Id:
x[Id_R] = 1.
if Id_L is not None:
x[Id_L] = 0.
x = npc.diag(x, wL_leg, labels=['wR*', 'wR'])
return x, separate_Id
class DMRGEngine(Sweep):
"""DMRG base class with common methods for the TwoSiteDMRG and SingleSiteDMRG.
This engine is implemented as a subclass of :class:`~tenpy.algorithms.mps_common.Sweep`.
It contains all methods that are generic between
:class:`SingleSiteDMRGEngine` and :class:`TwoSiteDMRGEngine`.
Use the latter two classes for actual DMRG runs.
A generic protocol for approaching a physics question using DMRG is given in :doc:`/intro/protocol`.
.. deprecated :: 0.5.0
Renamed parameter/attribute `DMRG_params` to :attr:`options`.
Options
-------
.. cfg:config :: DMRGEngine
:include: Sweep
Attributes
----------
EffectiveH : class type
Class for the effective Hamiltonian, i.e., a subclass of
:class:`~tenpy.algorithms.mps_common.EffectiveH`. Has a `length` class attribute which
specifies the number of sites updated at once (e.g., whether we do single-site vs. two-site
DMRG).
chi_list : dict | ``None``
See :cfg:option:`DMRGEngine.chi_list`
eff_H : :class:`~tenpy.algorithms.mps_common.EffectiveH`
Effective two-site Hamiltonian.
mixer : :class:`Mixer` | ``None``
If ``None``, no mixer is used (anymore), otherwise the mixer instance.
shelve : bool
If a simulation runs out of time (`time.time() - start_time > max_seconds`), the run will
terminate with `shelve = True`.
sweeps : int
The number of sweeps already performed. (Useful for re-start).
time0 : float
Time marker for the start of the run.
update_stats : dict
A dictionary with detailed statistics of the convergence at local update-level.
For each key in the following table, the dictionary contains a list where one value is
added each time :meth:`DMRGEngine.update_bond` is called.
=========== ===================================================================
key description
=========== ===================================================================
i0 An update was performed on sites ``i0, i0+1``.
----------- -------------------------------------------------------------------
age The number of physical sites involved in the simulation.
----------- -------------------------------------------------------------------
E_total The total energy before truncation.
----------- -------------------------------------------------------------------
N_lanczos Dimension of the Krylov space used in the lanczos diagonalization.
----------- -------------------------------------------------------------------
time Wallclock time evolved since :attr:`time0` (in seconds).
----------- -------------------------------------------------------------------
ov_change ``1. - abs(<theta_guess|theta_diag>)``, where ``|theta_guess>`` is
the initial guess for the wave function and ``|theta_diag>`` is the
*untruncated* wave function returned by :meth:`diag`.
=========== ===================================================================
sweep_stats : dict
A dictionary with detailed statistics at the sweep level.
For each key in the following table, the dictionary contains a list where one value is
added each time :meth:`Engine.sweep` is called (with ``optimize=True``).
============= ===================================================================
key description
============= ===================================================================
sweep Number of sweeps (excluding environment sweeps) performed so far.
------------- -------------------------------------------------------------------
N_updates Number of updates (including environment sweeps) performed so far.
------------- -------------------------------------------------------------------
E The energy *before* truncation (as calculated by Lanczos).
------------- -------------------------------------------------------------------
S Maximum entanglement entropy.
------------- -------------------------------------------------------------------
time Wallclock time evolved since :attr:`time0` (in seconds).
------------- -------------------------------------------------------------------
max_trunc_err The maximum truncation error in the last sweep
------------- -------------------------------------------------------------------
max_E_trunc Maximum change or Energy due to truncation in the last sweep.
------------- -------------------------------------------------------------------
max_chi Maximum bond dimension used.
------------- -------------------------------------------------------------------
norm_err Error of canonical form ``np.linalg.norm(psi.norm_test())``.
============= ===================================================================
"""
EffectiveH = None
DefaultMixer = None
def __init__(self, psi, model, options, **kwargs):
options = asConfig(options, self.__class__.__name__)
self.mixer = None
if isinstance(self, TwoSiteDMRGEngine):
self.DefaultMixer = DensityMatrixMixer
else:
self.DefaultMixer = SingleSiteMixer
super().__init__(psi, model, options, **kwargs)
@property
def DMRG_params(self):
warnings.warn("renamed self.DMRG_params -> self.options", FutureWarning, stacklevel=2)
return self.options
def run(self):
"""Run the DMRG simulation to find the ground state.
Returns
-------
E : float
The energy of the resulting ground state MPS.
psi : :class:`~tenpy.networks.mps.MPS`
The MPS representing the ground state after the simluation,
i.e. just a reference to :attr:`psi`.
Options
-------
.. cfg:configoptions :: DMRGEngine
diag_method : str
Method to be used for diagonalzation, default ``'default'``.
For possible arguments see :meth:`DMRGEngine.diag`.
E_tol_to_trunc : float
It's reasonable to choose the Lanczos convergence criteria
``'E_tol'`` not many magnitudes lower than the current
truncation error. Therefore, if `E_tol_to_trunc` is not
``None``, we update `E_tol` of `lanczos_params` to
``max_E_trunc*E_tol_to_trunc``,
restricted to the interval [`E_tol_min`, `E_tol_max`],
where ``max_E_trunc`` is the maximal energy difference due to
truncation right after each Lanczos optimization during the
sweeps.
E_tol_max : float
See `E_tol_to_trunc`
E_tol_min : float
See `E_tol_to_trunc`
max_E_err : float
Convergence if the change of the energy in each step
satisfies ``-Delta E / max(|E|, 1) < max_E_err``. Note that
this is also satisfied if ``Delta E > 0``,
i.e., if the energy increases (due to truncation).
max_hours : float
If the DMRG took longer (measured in wall-clock time),
'shelve' the simulation, i.e. stop and return with the flag
``shelve=True``.
max_S_err : float
Convergence if the relative change of the entropy in each step
satisfies ``|Delta S|/S < max_S_err``
max_sweeps : int
Maximum number of sweeps to be performed.
min_sweeps : int
Minimum number of sweeps to be performed.
Defaults to 1.5*N_sweeps_check.
N_sweeps_check : int
Number of sweeps to perform between checking convergence
criteria and giving a status update.
norm_tol : float
After the DMRG run, update the environment with at most
`norm_tol_iter` sweeps until
``np.linalg.norm(psi.norm_err()) < norm_tol``.
norm_tol_iter : float
Perform at most `norm_tol_iter`*`update_env` sweeps to
converge the norm error below `norm_tol`.
If the state is not converged after that, call
:meth:`~tenpy.networks.mps.canonical_form` instead.
P_tol_to_trunc : float
It's reasonable to choose the Lanczos convergence criteria
``'P_tol'`` not many magnitudes lower than the current
truncation error. Therefore, if `P_tol_to_trunc` is not
``None``, we update `P_tol` of `lanczos_params` to
``max_trunc_err*P_tol_to_trunc``,
restricted to the interval [`P_tol_min`, `P_tol_max`],
where ``max_trunc_err`` is the maximal truncation error
(discarded weight of the Schmidt values) due to truncation
right after each Lanczos optimization during the sweeps.
P_tol_max : float
See `P_tol_to_trunc`
P_tol_min : float
See `P_tol_to_trunc`
update_env : int
Number of sweeps without bond optimizaiton to update the
environment for infinite boundary conditions,
performed every `N_sweeps_check` sweeps.
"""
options = self.options
start_time = self.time0
self.shelve = False
# parameters for lanczos
p_tol_to_trunc = options.get('P_tol_to_trunc', 0.05)
if p_tol_to_trunc is not None:
p_tol_min = max(1.e-30,
self.trunc_params.silent_get('svd_min', 0.)**2 * p_tol_to_trunc,
self.trunc_params.silent_get('trunc_cut', 0.)**2 * p_tol_to_trunc)
p_tol_min = options.get('P_tol_min', p_tol_min)
p_tol_max = options.get('P_tol_max', 1.e-4)
e_tol_to_trunc = options.get('E_tol_to_trunc', None)
if e_tol_to_trunc is not None:
e_tol_min = options.get('E_tol_min', 5.e-16)
e_tol_max = options.get('E_tol_max', 1.e-4)
# parameters for DMRG convergence criteria
N_sweeps_check = options.get('N_sweeps_check', 1 if self.psi.finite else 10)
min_sweeps = int(1.5 * N_sweeps_check)
if self.chi_list is not None:
min_sweeps = max(max(self.chi_list.keys()), min_sweeps)
min_sweeps = options.get('min_sweeps', min_sweeps)
max_sweeps = options.get('max_sweeps', 1000)
max_E_err = options.get('max_E_err', 1.e-8)
max_S_err = options.get('max_S_err', 1.e-5)
max_seconds = 3600 * options.get('max_hours', 24 * 365)
norm_tol = options.get('norm_tol', 1.e-5)
if not self.finite:
update_env = options.get('update_env', N_sweeps_check // 2)
norm_tol_iter = options.get('norm_tol_iter', 5)
E_old, S_old = np.nan, np.mean(self.psi.entanglement_entropy()) # initial dummy values
E, Delta_E, Delta_S = 1., 1., 1.
self.diag_method = options.get('diag_method', 'default')
self.mixer_activate()
# loop over sweeps
while True:
loop_start_time = time.time()
# check convergence criteria
if self.sweeps >= max_sweeps:
break
if (self.sweeps > min_sweeps and -Delta_E < max_E_err * max(abs(E), 1.)
and abs(Delta_S) < max_S_err):
if self.mixer is None:
break
else:
logger.info("Convergence criterium reached with enabled mixer. "
"Disable mixer and continue")
self.mixer = None
if loop_start_time - start_time > max_seconds:
self.shelve = True
logger.warn("DMRG: maximum time limit reached. Shelve simulation.")
break
# --------- the main work --------------
logger.info('Running sweep with optimization')
for i in range(N_sweeps_check - 1):
self.sweep(meas_E_trunc=False)
max_trunc_err = self.sweep(meas_E_trunc=True)
max_E_trunc = np.max(self.E_trunc_list)
# --------------------------------------
# update lancos_params depending on truncation error(s)
if p_tol_to_trunc is not None and max_trunc_err > p_tol_min:
P_tol = max(p_tol_min, min(p_tol_max, max_trunc_err * p_tol_to_trunc))
self.lanczos_params['P_tol'] = P_tol
logger.debug("set lanczos_params['P_tol'] = %.2e", P_tol)
if e_tol_to_trunc is not None and max_E_trunc > e_tol_min:
E_tol = max(e_tol_min, min(e_tol_max, max_E_trunc * e_tol_to_trunc))
self.lanczos_params['E_tol'] = E_tol
logger.debug("set lanczos_params['E_tol'] = %.2e", E_tol)
# update environment
if not self.finite:
self.environment_sweeps(update_env)
# update values for checking the convergence
try:
S = self.psi.entanglement_entropy()
max_S = max(S)
S = np.mean(S)
Delta_S = (S - S_old) / N_sweeps_check
except ValueError:
# with a mixer, psi._S can be 2D arrays s.t. entanglement_entropy() fails
S = np.nan
max_S = np.nan
Delta_S = 0.
S_old = S
if not self.finite: # iDMRG: need energy density
Es = self.update_stats['E_total']
age = self.update_stats['age']
delta = min(1 + 2 * self.env.L, len(age))
growth = (age[-1] - age[-delta])
E = (Es[-1] - Es[-delta]) / growth
else:
E = self.update_stats['E_total'][-1]
Delta_E = (E - E_old) / N_sweeps_check
E_old = E
norm_err = np.linalg.norm(self.psi.norm_test())
# update statistics
self.sweep_stats['sweep'].append(self.sweeps)
self.sweep_stats['N_updates'].append(len(self.update_stats['i0']))
self.sweep_stats['E'].append(E)
self.sweep_stats['S'].append(S)
self.sweep_stats['time'].append(time.time() - start_time)
self.sweep_stats['max_trunc_err'].append(max_trunc_err)
self.sweep_stats['max_E_trunc'].append(max_E_trunc)
self.sweep_stats['max_chi'].append(np.max(self.psi.chi))
self.sweep_stats['norm_err'].append(norm_err)
# status update
logger.info(
"checkpoint after sweep %(sweeps)d\n"
"energy=%(E).16f, max S=%(S).16f, age=%(age)d, norm_err=%(norm_err).1e\n"
"Current memory usage %(mem).1fMB, wall time: %(wall_time).1fs\n"
"Delta E = %(dE).4e, Delta S = %(dS).4e (per sweep)\n"
"max trunc_err = %(trunc_err).4e, max E_trunc = %(E_trunc).4e\n"
"chi: %(chi)s\n"
"%(sep)s", {
'sweeps': self.sweeps,
'E': E,
'S': max_S,
'age': self.update_stats['age'][-1],
'norm_err': norm_err,
'mem': memory_usage(),
'wall_time': time.time() - loop_start_time,
'dE': Delta_E,
'dS': Delta_S,
'trunc_err': max_trunc_err,
'E_trunc': max_E_trunc,
'chi': self.psi.chi if self.psi.L < 40 else max(self.psi.chi),
'sep': "=" * 80,
})
self.checkpoint.emit(self)
# clean up from mixer
self.mixer_cleanup()
# update environment until norm_tol is reached
if norm_tol is not None and norm_err > norm_tol:
logger.warn(
"final DMRG state not in canonical form up to "
"norm_tol=%.2e: norm_err=%.2e", norm_tol, norm_err)
if self.finite:
self.psi.canonical_form()
else:
for _ in range(norm_tol_iter):
self.environment_sweeps(update_env)
norm_err = np.linalg.norm(self.psi.norm_test())
if norm_err <= norm_tol:
break
else:
logger.warn(
"norm_err=%.2e still too high after environment_sweeps, "
"call psi.canonical_form()", norm_err)
self.psi.canonical_form()
logger.info("DMRG finished after %d sweeps, max chi=%d", self.sweeps, max(self.psi.chi))
return E, self.psi
def reset_stats(self, resume_data=None):
"""Reset the statistics, useful if you want to start a new sweep run.
.. cfg:configoptions :: DMRGEngine
chi_list : dict | None
A dictionary to gradually increase the `chi_max` parameter of
`trunc_params`. The key defines starting from which sweep
`chi_max` is set to the value, e.g. ``{0: 50, 20: 100}`` uses
``chi_max=50`` for the first 20 sweeps and ``chi_max=100``
afterwards. Overwrites `trunc_params['chi_list']``.
By default (``None``) this feature is disabled.
sweep_0 : int
The number of sweeps already performed. (Useful for re-start).
"""
super().reset_stats(resume_data)
self.update_stats = {
'i0': [],
'age': [],
'E_total': [],
'N_lanczos': [],
'time': [],
'err': [],
'E_trunc': [],
'ov_change': []
}
self.sweep_stats = {
'sweep': [],