PurificationTEBD¶
full name: tenpy.algorithms.purification.PurificationTEBD
parent module:
tenpy.algorithms.purification
type: class
Inheritance Diagram
Methods
|
Initialize self. |
|
see |
|
Disentangle theta before splitting with svd. |
Try global disentangling by determining the maximally entangled pairs of sites. |
|
Perform a sweep through the system and disentangle with |
|
|
Generalization of |
Return necessary data to resume a |
|
Resume a run that was interrupted. |
|
Run TEBD real time evolution by N_steps`*`dt. |
|
TEBD algorithm in imaginary time to find the ground state. |
|
Run imaginary time evolution to cool down to the given beta. |
|
Returns list of necessary steps for the suzuki trotter decomposition. |
|
Return time steps of U for the Suzuki Trotter decomposition of desired order. |
|
|
Evolve by |
|
Updates the B matrices on a given bond. |
|
Update a bond with a (possibly non-unitary) U_bond. |
|
Perform an update suitable for imaginary time evolution. |
|
Updates either even or odd bonds in unit cell. |
Class Attributes and Properties
|
|
For each bond the total number of iterations performed in any |
|
truncation error introduced on each non-trivial bond. |
|
|
-
class
tenpy.algorithms.purification.
PurificationTEBD
(psi, model, options)[source]¶ Bases:
tenpy.algorithms.tebd.TEBDEngine
Time evolving block decimation (TEBD) for purification MPS.
Deprecated since version 0.6.0: Renamed parameter/attribute TEBD_params to
options
.- Parameters
psi (
PurificationMPS
) – Initial state to be time evolved. Modified in place.model (
NearestNeighborModel
) – The model representing the Hamiltonian for which we want to find the ground state.options (dict) – Further optional parameters as described in the following table. See
run()
andrun_GS()
for more details.
Options
-
config
PurificationTEBD
¶ option summary delta_tau_list (from TEBDEngine) in PurificationTEBD.run_GS
A list of floats: the timesteps to be used. [...]
Method for disentangling. [...]
dt (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Minimal time step by which to evolve.
N_steps (from TEBDEngine) in PurificationTEBD.run_GS
Number of steps before measurement can be performed
order (from TEBDEngine) in PurificationTEBD.run_GS
Order of the Suzuki-Trotter decomposition.
start_time (from TimeEvolutionAlgorithm) in TimeEvolutionAlgorithm
Initial value for :attr:`evolved_time`.
start_trunc_err (from TEBDEngine) in TEBDEngine
Initial truncation error for :attr:`trunc_err`.
trunc_params (from Algorithm) in Algorithm
Truncation parameters as described in :cfg:config:`truncation`.
-
option
disentangle
: None | str¶ Method for disentangling. See
get_disentangler()
.
-
option
-
used_disentangler
¶ The disentangler to be used on the auxiliar indices. Chosen by
get_disentangler()
, called with the TEBD parameter'disentangle'
. Defaults to the trivial disentangler foroptions['disentangle']=None
.- Type
Disentangler
-
_disent_iterations
¶ Number of iterations performed on all bonds, including trivial bonds; lenght L.
- Type
1D ndarray
-
_guess_U_disent
¶ Same index strucuture as self._U: for each two-site U of the physical time evolution the disentangler from the last application. Initialized to identities.
- Type
list of list of npc.Array
-
run_imaginary
(beta)[source]¶ Run imaginary time evolution to cool down to the given beta.
Note that we don’t change the norm attribute of the MPS, i.e. normalization is preserved.
- Parameters
beta (float) – The inverse temperature beta = 1/T, by which we should cool down. We evolve to the closest multiple of
options['dt']
, see alsoevolved_time
.
-
property
disent_iterations
¶ For each bond the total number of iterations performed in any
Disentangler
.
-
update_bond
(i, U_bond)[source]¶ Updates the B matrices on a given bond.
Function that updates the B matrices, the bond matrix s between and the bond dimension chi for bond i. This would look something like:
| | | | ... - B1 - s - B2 - ... | | | | |-------------| | | U | | |-------------| | | |
- Parameters
- Returns
trunc_err – The error of the represented state which is introduced by the truncation during this update step.
- Return type
-
update_bond_imag
(i, U_bond)[source]¶ Update a bond with a (possibly non-unitary) U_bond.
Similar as
update_bond()
; but after the SVD just keep the A, S, B canonical form. In that way, one can sweep left or right without using old singular values, thus preserving the canonical form during imaginary time evolution.- Parameters
- Returns
trunc_err – The error of the represented state which is introduced by the truncation during this update step.
- Return type
-
disentangle
(theta)[source]¶ Disentangle theta before splitting with svd.
For the purification we write \(\rho_P = Tr_Q{|\psi_{P,Q}><\psi_{P,Q}|}\). Thus, we can actually apply any unitary to the auxiliar Q space of \(|\psi>\) without changing the result.
Note
We have to apply the same unitary to the ‘bra’ and ‘ket’ used for expectation values / correlation functions!
The behaviour of this function is set by
used_disentangler
, which in turn is obtained fromget_disentangler(options['disentangle'])
, seeget_disentangler()
for details on the syntax.- Parameters
theta (
Array
) – Wave function to disentangle, with legs'vL', 'vR', 'p0', 'p1', 'q0', 'q1'
.- Returns
theta_disentangled (
Array
) – Disentangled theta;npc.tensordot(U, theta, axes=[['q0*', 'q1*'], ['q0', 'q1']])
.U (
Array
) – The unitary used to disentangle theta, with labels'q0', 'q1', 'q0*', 'q1*'
. If no unitary was found/applied, it might also beNone
.
-
disentangle_global
(pair=None)[source]¶ Try global disentangling by determining the maximally entangled pairs of sites.
Caclulate the mutual information (in the auxiliar space) between two sites and determine where it is maximal. Disentangle these two sites with
disentangle()
-
disentangle_global_nsite
(n=2)[source]¶ Perform a sweep through the system and disentangle with
disentangle_n_site()
.- Parameters
n (int) – maximal number of sites to disentangle at once.
-
disentangle_n_site
(i, n, theta)[source]¶ Generalization of
disentangle()
to n sites.Simply group left and right n/2 physical legs, adjust labels, and apply
disentangle()
to disentangle the central bond. Recursively proceed to disentangle left and right parts afterwards. Scales (for even n) as \(O(\chi^3 d^n d^{n/2})\).
-
get_resume_data
()[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(resume_data)
An algorithm which doesn’t support this should override resume_run to raise an Error.
- Returns
resume_data – Dictionary with necessary data (apart from copies of psi, model, options) that allows to continue the simulation from where we are now.
- 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 behaviour implemented here is to just callrun()
.
-
run_GS
()[source]¶ TEBD algorithm in imaginary time to find the ground state.
Note
It is almost always more efficient (and hence advisable) to use DMRG. This algorithms can nonetheless be used quite well as a benchmark and for comparison.
-
option
TEBDEngine
.
delta_tau_list
: list¶ A list of floats: the timesteps to be used. Choosing a large timestep delta_tau introduces large (Trotter) errors, but a too small time step requires a lot of steps to reach
exp(-tau H) --> |psi0><psi0|
. Therefore, we start with fairly large time steps for a quick time evolution until convergence, and the gradually decrease the time step.
-
option
TEBDEngine
.
order
: int¶ Order of the Suzuki-Trotter decomposition.
-
option
TEBDEngine
.
N_steps
: int¶ Number of steps before measurement can be performed
-
option
-
static
suzuki_trotter_decomposition
(order, N_steps)[source]¶ Returns list of necessary steps for the suzuki trotter decomposition.
We split the Hamiltonian as \(H = H_{even} + H_{odd} = H[0] + H[1]\). The Suzuki-Trotter decomposition is an approximation \(\exp(t H) \approx prod_{(j, k) \in ST} \exp(d[j] t H[k]) + O(t^{order+1 })\).
- Parameters
order (
1, 2, 4, '4_opt'
) – The desired order of the Suzuki-Trotter decomposition. Order1
approximation is simply \(e^A a^B\). Order2
is the “leapfrog” e^{A/2} e^B e^{A/2}. Order4
is the fourth-order from [suzuki1991] (also referenced in [schollwoeck2011]), and'4_opt'
gives the optmized version of Equ. (30a) in [barthel2020].- Returns
ST_decomposition – Indices
j, k
of the time-stepsd = suzuki_trotter_time_step(order)
and the decomposition of H. They are chosen such that a subsequent application ofexp(d[j] t H[k])
to a given state|psi>
yields(exp(N_steps t H[k]) + O(N_steps t^{order+1}))|psi>
.- Return type
-
static
suzuki_trotter_time_steps
(order)[source]¶ Return time steps of U for the Suzuki Trotter decomposition of desired order.
See
suzuki_trotter_decomposition()
for details.- Parameters
order (int) – The desired order of the Suzuki-Trotter decomposition.
- Returns
time_steps – We need
U = exp(-i H_{even/odd} delta_t * dt)
for the dt returned in this list.- Return type
list of float
-
property
trunc_err_bonds
¶ truncation error introduced on each non-trivial bond.
-
update
(N_steps)[source]¶ Evolve by
N_steps * U_param['dt']
.- Parameters
N_steps (int) – The number of steps for which the whole lattice should be updated.
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type
-
update_imag
(N_steps)[source]¶ Perform an update suitable for imaginary time evolution.
Instead of the even/odd brick structure used for ordinary TEBD, we ‘sweep’ from left to right and right to left, similar as DMRG. Thanks to that, we are actually able to preserve the canonical form.
- Parameters
N_steps (int) – The number of steps for which the whole lattice should be updated.
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type
-
update_step
(U_idx_dt, odd)[source]¶ Updates either even or odd bonds in unit cell.
Depending on the choice of p, this function updates all even (
E
, odd=False,0) or odd (O
) (odd=True,1) bonds:| - B0 - B1 - B2 - B3 - B4 - B5 - B6 - | | | | | | | | | | |----| |----| |----| | | | E | | E | | E | | | |----| |----| |----| | |----| |----| |----| | | | O | | O | | O | | | |----| |----| |----| |
Note that finite boundary conditions are taken care of by having
Us[0] = None
.- Parameters
U_idx_dt (int) – Time step index in
self._U
, evolve withUs[i] = self.U[U_idx_dt][i]
at bond(i-1,i)
.odd (bool/int) – Indication of whether to update even (
odd=False,0
) or even (odd=True,1
) sites
- Returns
trunc_err – The error of the represented state which is introduced due to the truncation during this sequence of update steps.
- Return type