EngineCombine¶
full name: tenpy.algorithms.dmrg.EngineCombine
parent module:
tenpy.algorithms.dmrg
type: class
-
class
tenpy.algorithms.dmrg.
EngineCombine
(psi, model, DMRG_params)[source]¶ Bases:
tenpy.algorithms.dmrg.TwoSiteDMRGEngine
Engine which combines legs into pipes as far as possible.
This engine combines the virtual and physical leg for the left site and right site into pipes. This reduces the overhead of calculating charge combinations in the contractions, but one
matvec()
is formally more expensive, \(O(2 d^3 \chi^3 D)\).Methods
diag
(self, theta_guess)Diagonalize the effective Hamiltonian represented by self.
environment_sweeps
(self, N_sweeps)Perform N_sweeps sweeps without optimization to update the environment.
get_sweep_schedule
(self)Define the schedule of the sweep.
init_env
(self[, model])(Re-)initialize the environment.
mixed_svd
(self, theta)Get (truncated) B from the new theta (as returned by diag).
mixer_activate
(self)Set self.mixer to the class specified by engine_params[‘mixer’].
mixer_cleanup
(self)Cleanup the effects of a mixer.
plot_sweep_stats
(self[, axes, xaxis, yaxis, …])Plot
sweep_stats
to display the convergence with the sweeps.plot_update_stats
(self, axes[, xaxis, …])Plot
update_stats
to display the convergence during the sweeps.post_update_local
(self, update_data[, …])Perform post-update actions.
prepare_svd
(self, theta)Transform theta into matrix for svd.
prepare_update
(self)Prepare self to represent the effective Hamiltonian on sites
(i0, i0+1)
.reset_stats
(self)Reset the statistics, useful if you want to start a new sweep run.
run
(self)Run the DMRG simulation to find the ground state.
set_B
(self, U, S, VH)Update the MPS with the
U, S, VH
returned by self.mixed_svd.sweep
(self[, optimize, meas_E_trunc])One ‘sweep’ of a sweeper algorithm.
update_LP
(self, U)Update left part of the environment.
update_RP
(self, VH)Update right part of the environment.
update_local
(self, theta[, optimize, …])Perform bond-update on the sites
(i0, i0+1)
.-
diag
(self, theta_guess)¶ Diagonalize the effective Hamiltonian represented by self.
The method used depends on the DMRG parameter diag_method.
diag_method
Function, comment
‘lanczos’
lanczos()
Default, the Lanczos implementation of TeNPy‘arpack’
lanczos_arpack()
Based onscipy.linalg.sparse.eigsh()
. Slower than ‘lanczos’, since it needs to convert the npc arrays to numpy arrays during each matvec, and possibly does many more iterations.‘ED_block’
full_diag_effH()
Contract the effective Hamiltonian to a (large!) matrix and diagonalize the block in the charge sector of the initial state. Preserves the charge sector of the explicitly conserved charges. However, if you don’t preserve a charge explicitly, it can break it. For example if you use aSpinChain({'conserve': 'parity'})
, it could change the total “Sz”, but not the parity of ‘Sz’.‘ED_all’
full_diag_effH()
Contract the effective Hamiltonian to a (large!) matrix and diagonalize it completely. Allows to change the charge sector even for explicitly conserved charges. For example if you use aSpinChain({'conserve': 'Sz'})
, it can change the total “Sz”.- Parameters
- theta_guess
Array
Initial guess for the ground state of the effective Hamiltonian.
- theta_guess
- Returns
- E0float
Energy of the found ground state.
- theta
Array
Ground state of the effective Hamiltonian.
- Nint
Number of Lanczos iterations used.
-1
if unknown.- ov_changefloat
Change in the wave function
1. - abs(<theta_guess|theta_diag>)
-
environment_sweeps
(self, N_sweeps)¶ Perform N_sweeps sweeps without optimization to update the environment.
- Parameters
- N_sweepsint
Number of sweeps to run without optimization
-
get_sweep_schedule
(self)¶ Define the schedule of the sweep.
One ‘sweep’ is a full sequence from the leftmost site to the right and back. Only those LP and RP that can be used later should be updated.
- Returns
- scheduleiterable of (int, bool, (bool, bool))
Schedule for the sweep. Each entry is
(i0, move_right, (update_LP, update_RP))
, where i0 is the leftmost of theself.EffectiveH.length
sites to be updated inupdate_local()
, move_right indicates whether the next i0 in the schedule is rigth (True) of the current one, and update_LP, update_RP indicate whether it is necessary to update the LP and RP. The latter are chosen such that the environment is growing for infinite systems, but we only keep the minimal number of environment tensors in memory.
-
init_env
(self, model=None)¶ (Re-)initialize the environment.
This function is useful to (re-)start a Sweep with a slightly different model or different (engine) parameters. Note that we assume that we still have the same psi. Calls
reset_stats()
.- Parameters
- model
MPOModel
The model representing the Hamiltonian for which we want to find the ground state. If
None
, keep the model used before.
- model
- Raises
- ValueError
If the engine is re-initialized with a new model, which legs are incompatible with those of hte old model.
-
mixed_svd
(self, theta)¶ Get (truncated) B from the new theta (as returned by diag).
The goal is to split theta and truncate it:
| -- theta -- ==> -- U -- S -- VH - | | | | |
Without a mixer, this is done by a simple svd and truncation of Schmidt values.
With a mixer, the state is perturbed before the SVD. The details of the perturbation are defined by the
Mixer
class.Note that the returned S is a general (not diagonal) matrix, with labels
'vL', 'vR'
.- Parameters
- theta
Array
The optimized wave function, prepared for svd.
- theta
- Returns
- U
Array
Left-canonical part of theta. Labels
'(vL.p0)', 'vR'
.- S1D ndarray | 2D
Array
Without mixer just the singluar values of the array; with mixer it might be a general matrix with labels
'vL', 'vR'
; see comment above.- VH
Array
Right-canonical part of theta. Labels
'vL', '(p1.vR)'
.- err
TruncationError
The truncation error introduced.
- U
-
mixer_activate
(self)¶ Set self.mixer to the class specified by engine_params[‘mixer’].
-
mixer_cleanup
(self)¶ Cleanup the effects of a mixer.
A
sweep()
with an enabledMixer
leaves the MPS psi with 2D arrays in S. To recover the originial form, this function simply performs one sweep with disabled mixer.
-
plot_sweep_stats
(self, axes=None, xaxis='time', yaxis='E', y_exact=None, **kwargs)¶ Plot
sweep_stats
to display the convergence with the sweeps.- Parameters
- axes
matplotlib.axes.Axes
The axes to plot into. Defaults to
matplotlib.pyplot.gca()
- xaxis, yaxiskey of
sweep_stats
Key of
sweep_stats
to be used for the x-axis and y-axis of the plots.- y_exactfloat
Exact value for the quantity on the y-axis for comparison. If given, plot
abs((y-y_exact)/y_exact)
on a log-scale yaxis.- **kwargs :
Further keyword arguments given to
axes.plot(...)
.
- axes
-
plot_update_stats
(self, axes, xaxis='time', yaxis='E', y_exact=None, **kwargs)¶ Plot
update_stats
to display the convergence during the sweeps.- Parameters
- axes
matplotlib.axes.Axes
The axes to plot into. Defaults to
matplotlib.pyplot.gca()
- xaxis
'N_updates' | 'sweep'
| keys ofupdate_stats
Key of
update_stats
to be used for the x-axis of the plots.'N_updates'
is just enumerating the number of bond updates, and'sweep'
corresponds to the sweep number (including environment sweeps).- yaxis
'E'
| keys ofupdate_stats
Key of
update_stats
to be used for the y-axisof the plots. For ‘E’, use the energy (per site for infinite systems).- y_exactfloat
Exact value for the quantity on the y-axis for comparison. If given, plot
abs((y-y_exact)/y_exact)
on a log-scale yaxis.- **kwargs :
Further keyword arguments given to
axes.plot(...)
.
- axes
-
post_update_local
(self, update_data, meas_E_trunc=False)¶ Perform post-update actions.
Compute truncation energy, remove LP/RP that are no longer needed and collect statistics.
- Parameters
- update_datadict
Data computed during the local update, as described in the following list.
- meas_E_truncbool, optional
Wheter to measure the energy after truncation.
-
prepare_svd
(self, theta)¶ Transform theta into matrix for svd.
-
prepare_update
(self)¶ Prepare self to represent the effective Hamiltonian on sites
(i0, i0+1)
.- Returns
- theta
Array
Current best guess for the ground state, which is to be optimized. Labels
'vL', 'p0', 'vR', 'p1'
.
- theta
-
reset_stats
(self)¶ Reset the statistics, useful if you want to start a new sweep run.
-
run
(self)¶ Run the DMRG simulation to find the ground state.
- Returns
- Efloat
The energy of the resulting ground state MPS.
- psi
MPS
The MPS representing the ground state after the simluation, i.e. just a reference to
psi
.
-
set_B
(self, U, S, VH)¶ Update the MPS with the
U, S, VH
returned by self.mixed_svd.
-
sweep
(self, optimize=True, meas_E_trunc=False)¶ One ‘sweep’ of a sweeper algorithm.
Iteratate over the bond which is optimized, to the right and then back to the left to the starting point. If optimize=False, don’t actually diagonalize the effective hamiltonian, but only update the environment.
- Parameters
- optimizebool, optional
Whether we actually optimize to find the ground state of the effective Hamiltonian. (If False, just update the environments).
- meas_E_truncbool, optional
Whether to measure truncation energies.
- Returns
- max_trunc_errfloat
Maximal truncation error introduced.
- max_E_trunc
None
| float None
if meas_E_trunc is False, else the maximal change of the energy due to the truncation.
-
update_LP
(self, U)¶ Update left part of the environment.
We always update the environment at site i0 + 1: this environment then contains the site where we just performed a local update (when sweeping right).
- Parameters
- U
Array
The U as returned by the SVD, with combined legs, labels
'vL.p0', 'vR'
.
- U
-
update_RP
(self, VH)¶ Update right part of the environment.
We always update the environment at site i0: this environment then contains the site where we just performed a local update (when sweeping left).
- Parameters
- VH
Array
The VH as returned by SVD, with combined legs, labels
'vL', '(vR.p1)'
.
- VH
-
update_local
(self, theta, optimize=True, meas_E_trunc=False)¶ Perform bond-update on the sites
(i0, i0+1)
.- Parameters
- theta
Array
Initial guess for the ground state of the effective Hamiltonian.
- optimizebool
Wheter we actually optimize to find the ground state of the effective Hamiltonian. (If False, just update the environments).
- meas_E_truncbool
Wheter to measure the energy after truncation.
- theta
- Returns
- update_datadict
Data computed during the local update, as described in the following:
- E0float
Total energy, obtained before truncation (if
optimize=True
), or after truncation (ifoptimize=False
) (but neverNone
).- Nint
Dimension of the Krylov space used for optimization in the lanczos algorithm. 0 if
optimize=False
.- ageint
Current size of the DMRG simulation: number of physical sites involved into the contraction.
- U, VH:
Array
U and VH returned by
mixed_svd()
.- ov_change: float
Change in the wave function
1. - abs(<theta_guess|theta>)
induced bydiag()
, not including the truncation!
-