DMRGEngine
full name: tenpy.algorithms.dmrg.DMRGEngine
parent module:
tenpy.algorithms.dmrg
type: class
Inheritance Diagram
Methods
|
|
|
Diagonalize the effective Hamiltonian represented by self. |
|
Perform N_sweeps sweeps without optimization to update the environment. |
|
Gives an approximate prediction for the required memory usage. |
Remove no longer needed environments after an update. |
|
|
Return necessary data to resume a |
Define the schedule of the sweep. |
|
|
(Re-)initialize the environment. |
Determines if the algorithm is converged. |
|
Create new instance of self.EffectiveH at self.i0 and set it to self.eff_H. |
|
Set self.mixer to the class specified by options['mixer']. |
|
Cleanup the effects of a mixer. |
|
Deactivate the mixer. |
|
|
Plot |
|
Plot |
Perform any final steps or clean up after the main loop has terminated. |
|
|
Perform post-update actions. |
Perform preparations before |
|
Prepare self for calling |
|
|
Reset the statistics, useful if you want to start a new sweep run. |
Resume a run that was interrupted. |
|
Run the DMRG simulation to find the ground state. |
|
Perform a single iteration, consisting of |
|
|
Emits a status message to the logging system after an iteration. |
Determines if the main loop should be terminated. |
|
|
One 'sweep' of the algorithm. |
|
Initialize algorithm from another algorithm instance of a different class. |
|
Update the left and right environments after an update of the state. |
|
Perform site-update on the site |
Update the singular values at the boundaries of the segment. |
Class Attributes and Properties
|
|
|
|
|
|
The number of sites to be optimized at once. |
|
|
- class tenpy.algorithms.dmrg.DMRGEngine(psi, model, options, **kwargs)[source]
Bases:
IterativeSweeps
DMRG base class with common methods for the TwoSiteDMRG and SingleSiteDMRG.
This engine is implemented as a subclass of
Sweep
. It contains all methods that are generic betweenSingleSiteDMRGEngine
andTwoSiteDMRGEngine
. Use the latter two classes for actual DMRG runs.A generic protocol for approaching a physics question using DMRG is given in Protocol for using (i)DMRG.
Options
- config DMRGEngine
option summary By default (``None``) this feature is disabled. [...]
chi_list_reactivates_mixer (from Sweep) in IterativeSweeps.sweep
If True, the mixer is reset/reactivated each time the bond dimension growth [...]
Whether to combine legs into pipes. This combines the virtual and [...]
diag_method in DMRGEngine.diag
One of the following strings: [...]
E_tol_max in DMRGEngine.run_iteration
See `E_tol_to_trunc`
E_tol_min in DMRGEngine.run_iteration
See `E_tol_to_trunc`
E_tol_to_trunc in DMRGEngine.run_iteration
It's reasonable to choose the Lanczos convergence criteria [...]
lanczos_params (from Sweep) in Sweep
Lanczos parameters as described in :cfg:config:`KrylovBased`.
max_E_err in DMRGEngine.is_converged
Convergence if the change of the energy in each step [...]
max_hours (from IterativeSweeps) in DMRGEngine.stopping_criterion
If the DMRG took longer (measured in wall-clock time), [...]
max_N_for_ED in DMRGEngine.diag
Maximum matrix dimension of the effective hamiltonian [...]
max_N_sites_per_ring (from Algorithm) in Algorithm
Threshold for raising errors on too many sites per ring. Default ``18``. [...]
max_S_err in DMRGEngine.is_converged
Convergence if the relative change of the entropy in each step [...]
max_sweeps (from IterativeSweeps) in DMRGEngine.stopping_criterion
Maximum number of sweeps to perform.
max_trunc_err (from IterativeSweeps) in IterativeSweeps
Threshold for raising errors on too large truncation errors. Default ``0.00 [...]
min_sweeps (from IterativeSweeps) in DMRGEngine.stopping_criterion
Minimum number of sweeps to perform.
Specifies which :class:`Mixer` to use, if any. [...]
mixer_params (from Sweep) in DMRGEngine.mixer_activate
Mixer parameters as described in :cfg:config:`Mixer`.
N_sweeps_check in DMRGEngine.run_iteration
Number of sweeps to perform between checking convergence [...]
norm_tol in DMRGEngine.post_run_cleanup
After the DMRG run, update the environment with at most [...]
norm_tol_final in DMRGEngine.post_run_cleanup
After performing `norm_tol_iter`*`update_env` sweeps, if [...]
norm_tol_iter in DMRGEngine.post_run_cleanup
Perform at most `norm_tol_iter`*`update_env` sweeps to [...]
P_tol_max in DMRGEngine.run_iteration
See `P_tol_to_trunc`
P_tol_min in DMRGEngine.run_iteration
See `P_tol_to_trunc`
P_tol_to_trunc in DMRGEngine.run_iteration
It's reasonable to choose the Lanczos convergence criteria [...]
Number of sweeps to be performed without optimization to update the environment.
trunc_params (from Algorithm) in Algorithm
Truncation parameters as described in :cfg:config:`truncation`.
update_env in DMRGEngine.run_iteration
Number of sweeps without bond optimization to update the [...]
- update_stats
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
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
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 bydiag()
.- Type:
- sweep_stats
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
Engine.sweep()
is called (withoptimize=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).
Delta_E
The change in E (above) since the last iteration.
S
Mean entanglement entropy (over bonds).
Delta_S
The change in S (above) since the last iteration.
max_S
Max entanglement entropy (over bonds).
time
Wallclock time evolved since
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())
.- Type:
- _entropy_approx
While the mixer is on, the S stored in the MPS is a non-diagonal 2D array. To check convergence, we use the approximate singular values based on which we truncated instead to calculate the entanglement entropy and store it inside this list.
- Type:
list of {None, 1D array}
- pre_run_initialize()[source]
Perform preparations before
run_iteration()
is iterated.- Returns:
The object to be returned by
run()
in case of immediate convergence, i.e. if no iterations are performed.- Return type:
result
- run_iteration()[source]
Perform a single iteration, consisting of
N_sweeps_check
sweeps.Options
- option DMRGEngine.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 notNone
, we update E_tol of lanczos_params tomax_E_trunc*E_tol_to_trunc
, restricted to the interval [E_tol_min, E_tol_max], wheremax_E_trunc
is the maximal energy difference due to truncation right after each Lanczos optimization during the sweeps.
- option DMRGEngine.E_tol_max: float
See E_tol_to_trunc
- option DMRGEngine.E_tol_min: float
See E_tol_to_trunc
- option DMRGEngine.N_sweeps_check: int
Number of sweeps to perform between checking convergence criteria and giving a status update.
- option DMRGEngine.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 notNone
, we update P_tol of lanczos_params tomax_trunc_err*P_tol_to_trunc
, restricted to the interval [P_tol_min, P_tol_max], wheremax_trunc_err
is the maximal truncation error (discarded weight of the Schmidt values) due to truncation right after each Lanczos optimization during the sweeps.
- option DMRGEngine.P_tol_max: float
See P_tol_to_trunc
- option DMRGEngine.P_tol_min: float
See P_tol_to_trunc
- option DMRGEngine.update_env: int
Number of sweeps without bond optimization to update the environment for infinite boundary conditions, performed every N_sweeps_check sweeps.
- Returns:
E (float) – The energy of the current ground state approximation.
psi (
MPS
) – The current ground state approximation, i.e. just a reference topsi
.
- status_update(iteration_start_time: float)[source]
Emits a status message to the logging system after an iteration.
- Parameters:
iteration_start_time (float) – The
time.time()
at the start of the last iteration
- is_converged()[source]
Determines if the algorithm is converged.
Does not cover any other reasons to abort, such as reaching a time limit. Such checks are covered by
stopping_condition()
.Options
- option DMRGEngine.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 might be satisfied even ifDelta E > 0
, i.e., if the energy increases (due to truncation).
- option DMRGEngine.max_S_err: float
Convergence if the relative change of the entropy in each step satisfies
|Delta S|/S < max_S_err
- post_run_cleanup()[source]
Perform any final steps or clean up after the main loop has terminated.
Options
- option DMRGEngine.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
.
- option DMRGEngine.norm_tol_iter: float
Perform at most norm_tol_iter`*`update_env sweeps to converge the norm error below norm_tol.
- option DMRGEngine.norm_tol_final: float
After performing norm_tol_iter`*`update_env sweeps, if
np.linalg.norm(psi.norm_err()) < norm_tol_final
, callcanonical_form()
to canonicalize instead. This tolerance should be stricter than norm_tol to ensure canonical form even if DMRG cannot fully converge.
- run()[source]
Run the DMRG simulation to find the ground state.
- Returns:
E (float) – The energy of the resulting ground state MPS.
psi (
MPS
) – The MPS representing the ground state after the simulation, i.e. just a reference topsi
.
- reset_stats(resume_data=None)[source]
Reset the statistics, useful if you want to start a new sweep run.
- sweep(optimize=True, meas_E_trunc=False)[source]
One ‘sweep’ of the algorithm.
Thin wrapper around
tenpy.algorithms.mps_common.Sweep.sweep()
with one additional parameter meas_E_trunc specifying whether to measure truncation energies.
- update_local(theta, optimize=True)[source]
Perform site-update on the site
i0
.- Parameters:
- Returns:
update_data – 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!
- Return type:
- post_update_local(E0, age, N, ov_change, err, **update_data)[source]
Perform post-update actions.
Compute truncation energy and collect statistics.
- Parameters:
**update_data (dict) – What was returned by
update_local()
.
- update_segment_boundaries()[source]
Update the singular values at the boundaries of the segment.
This method is called at the end of
post_update_local()
for ‘segment’ boundary MPS. It just updates the singular values on the very left/right end of the MPS segment.
- diag(theta_guess)[source]
Diagonalize the effective Hamiltonian represented by self.
- option DMRGEngine.max_N_for_ED: int
Maximum matrix dimension of the effective hamiltonian up to which the
'default'
diag_method uses ED instead of Lanczos.
- option DMRGEngine.diag_method: str
One of the following strings:
- ‘default’
Same as
'lanczos'
for large bond dimensions, but if the total dimension of the effective Hamiltonian does not exceed the DMRG parameter'max_N_for_ED'
it uses'ED_block'
.- ‘lanczos’
lanczos()
Default, the Lanczos implementation in 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.- Returns:
E0 (float) – Energy of the found ground state.
theta (
Array
) – Ground state of the effective Hamiltonian.N (int) – Number of Lanczos iterations used.
-1
if unknown.ov_change (float) – Change in the wave function
1. - abs(<theta_guess|theta_diag>)
- plot_update_stats(axes, xaxis='time', yaxis='E', y_exact=None, **kwargs)[source]
Plot
update_stats
to display the convergence during the sweeps.- Parameters:
axes (
matplotlib.axes.Axes
) – The axes to plot into. Defaults tomatplotlib.pyplot.gca()
xaxis (
'N_updates' | 'sweep'
| keys ofupdate_stats
) – Key ofupdate_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 ofupdate_stats
to be used for the y-axis of the plots. For ‘E’, use the energy (per site for infinite systems).y_exact (float) – 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(...)
.
- plot_sweep_stats(axes=None, xaxis='time', yaxis='E', y_exact=None, **kwargs)[source]
Plot
sweep_stats
to display the convergence with the sweeps.- Parameters:
axes (
matplotlib.axes.Axes
) – The axes to plot into. Defaults tomatplotlib.pyplot.gca()
xaxis (key of
sweep_stats
) – Key ofsweep_stats
to be used for the x-axis and y-axis of the plots.yaxis (key of
sweep_stats
) – Key ofsweep_stats
to be used for the x-axis and y-axis of the plots.y_exact (float) – 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(...)
.
- environment_sweeps(N_sweeps)[source]
Perform N_sweeps sweeps without optimization to update the environment.
- Parameters:
N_sweeps (int) – Number of sweeps to run without optimization
- estimate_RAM(mem_saving_factor=None)[source]
Gives an approximate prediction for the required memory usage.
This calculation is based on the requested bond dimension, the local Hilbert space dimension, the number of sites, and the boundary conditions.
- Parameters:
mem_saving_factor (float) – Represents the amount of RAM saved due to conservation laws. By default, it is ‘None’ and is extracted from the model automatically. However, this is only possible in a few cases and needs to be estimated in most cases. This is due to the fact that it is dependent on the model parameters. If one has a better estimate, one can pass the value directly. This value can be extracted by building the initial state psi (usually by performing DMRG) and then calling
print(psi.get_B(0).sparse_stats())
TeNPy will automatically print the fraction of nonzero entries in the first line, for example,6 of 16 entries (=0.375) nonzero
. This fraction corresponds to the mem_saving_factor; in our example, it is 0.375.- Returns:
usage – Required RAM in MB.
- Return type:
See also
tenpy.simulations.simulation.estimate_simulation_RAM
global function calling this.
- free_no_longer_needed_envs()[source]
Remove no longer needed environments after an update.
This allows to minimize the number of environments to be kept. For large MPO bond dimensions, these environments are by far the biggest part in memory, so this is a valuable optimization to reduce memory requirements.
- get_resume_data(sequential_simulations=False)[source]
Return necessary data to resume a
run()
interrupted at a checkpoint.At a
checkpoint
, you can savepsi
,model
andoptions
along with the data returned by this function. When the simulation aborts, you can resume it using this saved data with:eng = AlgorithmClass(psi, model, options, resume_data=resume_data) eng.resume_run()
An algorithm which doesn’t support this should override resume_run to raise an Error.
- Parameters:
sequential_simulations (bool) – If True, return only the data for re-initializing a sequential simulation run, where we “adiabatically” follow the evolution of a ground state (for variational algorithms), or do series of quenches (for time evolution algorithms); see
run_seq_simulations()
.- Returns:
resume_data – Dictionary with necessary data (apart from copies of psi, model, options) that allows to continue the algorithm run from where we are now. It might contain an explicit copy of psi.
- Return type:
- get_sweep_schedule()[source]
Define the schedule of the sweep.
One ‘sweep’ is a full sequence from the leftmost site to the right and back.
- Returns:
schedule – 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 right (True), left (False) or equal (None) of the current one, and update_LP, update_RP indicate whether it is necessary to update the LP and RP of the environments.- Return type:
- init_env(model=None, resume_data=None, orthogonal_to=None)[source]
(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. IfNone
, keep the model used before.resume_data (None | dict) – Given when resuming a simulation, as returned by
get_resume_data()
. Can contain another dict under the key init_env_data; the contents of init_env_data get passed as keyword arguments to the environment initialization.orthogonal_to (None | list of
MPS
| list of dict) – List of other matrix product states to orthogonalize against. Instead of just the state, you can specify a dict with the state as ket and further keyword arguments for initializing theMPSEnvironment
; thepsi
to be optimized is used as bra. Works only for finite or segment MPS; for infinite MPS it must be None. This can be used to find (a few) excited states as follows. First, run DMRG to find the ground state, and then run DMRG again while orthogonalizing against the ground state, which yields the first excited state (in the same symmetry sector), and so on. Note thatresume_data['orthogonal_to']
takes precedence over the argument.
Options
- option Sweep.start_env: int
Number of sweeps to be performed without optimization to update the environment.
- Raises:
ValueError – If the engine is re-initialized with a new model, which legs are incompatible with those of hte old model.
- mixer_activate()[source]
Set self.mixer to the class specified by options[‘mixer’].
- option Sweep.mixer: str | class | bool | None
Specifies which
Mixer
to use, if any. A string stands for one of the mixers defined in this module. A class is assumed to have the same interface asMixer
and is used to instantiate themixer
.None
uses no mixer.True
uses the mixer specified by theDefaultMixer
class attribute. The default depends on the subclass ofSweep
.
See also
- mixer_cleanup()[source]
Cleanup the effects of a mixer.
A
sweep()
with an enabledMixer
leaves the MPS psi with 2D arrays in S. This method recovers the original form by performing SVDs of the S and updating the MPS tensors accordingly.
- mixer_deactivate()[source]
Deactivate the mixer.
Set
self.mixer=None
and revert any other effects ofmixer_activate()
.
- property n_optimize
The number of sites to be optimized at once.
Indirectly set by the class attribute
EffectiveH
and it’s length. For example,TwoSiteDMRGEngine
uses theTwoSiteH
and hence hasn_optimize=2
, while theSingleSiteDMRGEngine
hasn_optimize=1
.
- prepare_update_local()[source]
Prepare self for calling
update_local()
.- Returns:
theta – Current best guess for the ground state, which is to be optimized. Labels are
'vL', 'p0', 'p1', 'vR'
, or combined versions of it (if self.combine). For single-site DMRG, the'p1'
label is missing.- Return type:
- resume_run()[source]
Resume a run that was interrupted.
In case we saved an intermediate result at a
checkpoint
, this function allows to resume therun()
of the algorithm (after re-initialization with the resume_data). Since most algorithms just have a while loop with break conditions, the default behavior implemented here is to just callrun()
.
- stopping_criterion(iteration_start_time: float) bool [source]
Determines if the main loop should be terminated.
- Parameters:
iteration_start_time (float) – The
time.time()
at the start of the last iteration
Options
- option IterativeSweeps.min_sweeps: int
Minimum number of sweeps to perform.
- option IterativeSweeps.max_sweeps: int
Maximum number of sweeps to perform.
- option IterativeSweeps.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
.
- classmethod switch_engine(other_engine, *, options=None, **kwargs)[source]
Initialize algorithm from another algorithm instance of a different class.
You can initialize one engine from another, not too different subclasses. Internally, this function calls
get_resume_data()
to extract data from the other_engine and then initializes the new class.Note that it transfers the data without making copies in most case; even the options! Thus, when you call run() on one of the two algorithm instances, it will modify the state, environment, etc. in the other. We recommend to make the switch as
engine = OtherSubClass.switch_engine(engine)
directly replacing the reference.- Parameters:
cls (class) – Subclass of
Algorithm
to be initialized.other_engine (
Algorithm
) – The engine from which data should be transferred. Another, but not too different algorithm subclass-class; e.g. you can switch from theTwoSiteDMRGEngine
to theOneSiteDMRGEngine
.options (None | dict-like) – If not None, these options are used for the new initialization. If None, take the options from the other_engine.
**kwargs – Further keyword arguments for class initialization. If not defined, resume_data is collected with
get_resume_data()
.
- update_env(**update_data)[source]
Update the left and right environments after an update of the state.
- Parameters:
**update_data – Whatever is returned by
update_local()
.