Engine¶
full name: tenpy.algorithms.tdvp.Engine
parent module:
tenpy.algorithms.tdvp
type: class
-
class
tenpy.algorithms.tdvp.
Engine
(psi, model, TDVP_params, environment=None)[source]¶ Bases:
object
Time dependant variational principle ‘Engine’.
You can call
run_one_site()
for single-site TDVP, orrun_two_sites()
for two-site TDVP.- Parameters
- psi
MPS
Initial state to be time evolved. Modified in place.
- model
MPOModel
The model representing the Hamiltonian for which we want to find the ground state.
- TDVP_paramsdict
Further optional parameters as described in the following table. Use
verbose>0
to print the used parameters during runtime.key
type
description
start_time
float
Initial value for
evolved_time
dt
float
Time step of the Trotter error
trunc_params
dict
Truncation parameters as described in
truncate()
- environment:class:’~tenpy.networks.mpo.MPOEnvironment` | None
Initial environment. If
None
(default), it will be calculated at the beginning.
- psi
- Attributes
- verboseint
Level of verbosity (i.e. how much status information to print); higher=more output.
- evolved_timefloat | complex
Indicating how long psi has been evolved,
psi = exp(-i * evolved_time * H) psi(t=0)
.- psi
MPS
The MPS, time evolved in-place.
- TDVP_params: dict
Optional parameters, see
run()
andrun_GS()
for more details.- environment
MPOEnvironment
The environment, storing the LP and RP to avoid recalculations.
Methods
run_one_site
(self[, N_steps])Run the TDVP algorithm with the one site algorithm.
run_two_sites
(self[, N_steps])Run the TDVP algorithm with two sites update.
set_anonymous_svd
(self, U, new_label)Relabel the svd.
sweep_left_right
(self)Performs the sweep left->right of the second order TDVP scheme with one site update.
sweep_left_right_two
(self)Performs the sweep left->right of the second order TDVP scheme with two sites update.
sweep_right_left
(self)Performs the sweep right->left of the second order TDVP scheme with one site update.
sweep_right_left_two
(self)Performs the sweep left->right of the second order TDVP scheme with two sites update.
theta_svd_left_right
(self, theta)Performs the SVD from left to right.
theta_svd_right_left
(self, theta)Performs the SVD from right to left.
update_s_h0
(self, s, H, dt)Update with the zero site Hamiltonian (update of the singular value)
update_theta_h1
(self, Lp, Rp, theta, W, dt)Update with the one site Hamiltonian.
update_theta_h2
(self, Lp, Rp, theta, W0, W1, dt)Update with the two sites Hamiltonian.
-
run_one_site
(self, N_steps=None)[source]¶ Run the TDVP algorithm with the one site algorithm.
Warning
Be aware that the bond dimension will not increase!
- Parameters
- N_stepsinteger. Number of steps
-
run_two_sites
(self, N_steps=None)[source]¶ Run the TDVP algorithm with two sites update.
The bond dimension will increase. Truncation happens at every step of the sweep, according to the parameters set in trunc_params.
- Parameters
- N_stepsinteger. Number of steps
-
sweep_left_right
(self)[source]¶ Performs the sweep left->right of the second order TDVP scheme with one site update.
Evolve from 0.5*dt.
-
sweep_left_right_two
(self)[source]¶ Performs the sweep left->right of the second order TDVP scheme with two sites update.
Evolve from 0.5*dt
-
sweep_right_left
(self)[source]¶ Performs the sweep right->left of the second order TDVP scheme with one site update.
Evolve from 0.5*dt
-
sweep_right_left_two
(self)[source]¶ Performs the sweep left->right of the second order TDVP scheme with two sites update.
Evolve from 0.5*dt
-
update_theta_h2
(self, Lp, Rp, theta, W0, W1, dt)[source]¶ Update with the two sites Hamiltonian.
- Parameters
- Lp
tenpy.linalg.np_conserved.Array
tensor representing the left environment
- Rp
tenpy.linalg.np_conserved.Array
tensor representing the right environment
- theta
tenpy.linalg.np_conserved.Array
the theta tensor which needs to be updated
- W
tenpy.linalg.np_conserved.Array
MPO which is applied to the ‘p0’ leg of theta
- W1
tenpy.linalg.np_conserved.Array
MPO which is applied to the ‘p1’ leg of theta
- Lp
-
theta_svd_left_right
(self, theta)[source]¶ Performs the SVD from left to right.
- Parameters
- theta: :class:`tenpy.linalg.np_conserved.Array`
the theta tensor on which the SVD is applied
-
set_anonymous_svd
(self, U, new_label)[source]¶ Relabel the svd.
- Parameters
- U
tenpy.linalg.np_conserved.Array
the tensor which lacks a leg_label
- U
-
theta_svd_right_left
(self, theta)[source]¶ Performs the SVD from right to left.
- Parameters
- theta
tenpy.linalg.np_conserved.Array
, The theta tensor on which the SVD is applied
- theta
-
update_s_h0
(self, s, H, dt)[source]¶ Update with the zero site Hamiltonian (update of the singular value)
- Parameters
- s
tenpy.linalg.np_conserved.Array
representing the singular value matrix which is updated
- HH0_mixed
zero site Hamiltonian that we need to apply on the singular value matrix
- dtcomplex number
time step of the evolution
- s