TimeDependentTwoSiteTDVP

  • full name: tenpy.algorithms.tdvp.TimeDependentTwoSiteTDVP

  • parent module: tenpy.algorithms.tdvp

  • type: class

Inheritance Diagram

Inheritance diagram of tenpy.algorithms.tdvp.TimeDependentTwoSiteTDVP

Methods

TimeDependentTwoSiteTDVP.__init__(psi, ...)

TimeDependentTwoSiteTDVP.environment_sweeps(...)

Perform N_sweeps sweeps without optimization to update the environment.

TimeDependentTwoSiteTDVP.evolve(N_steps, dt)

Evolve by N_steps * dt.

TimeDependentTwoSiteTDVP.evolve_step(dt)

TimeDependentTwoSiteTDVP.free_no_longer_needed_envs()

Remove no longer needed environments after an update.

TimeDependentTwoSiteTDVP.get_resume_data([...])

Return necessary data to resume a run() interrupted at a checkpoint.

TimeDependentTwoSiteTDVP.get_sweep_schedule()

Slightly different sweep schedule than DMRG

TimeDependentTwoSiteTDVP.init_env([model, ...])

(Re-)initialize the environment.

TimeDependentTwoSiteTDVP.make_eff_H()

Create new instance of self.EffectiveH at self.i0 and set it to self.eff_H.

TimeDependentTwoSiteTDVP.one_site_update(i, dt)

TimeDependentTwoSiteTDVP.post_update_local(...)

Algorithm-specific actions to be taken after local update.

TimeDependentTwoSiteTDVP.prepare_evolve(dt)

Do nothing.

TimeDependentTwoSiteTDVP.prepare_update_local()

Prepare self for calling update_local().

TimeDependentTwoSiteTDVP.reinit_model()

Re-initialize a new model at current evolved_time.

TimeDependentTwoSiteTDVP.reset_stats([...])

Reset the statistics.

TimeDependentTwoSiteTDVP.resume_run()

Resume a run that was interrupted.

TimeDependentTwoSiteTDVP.run()

Perform a (real-)time evolution of psi by N_steps * dt.

TimeDependentTwoSiteTDVP.run_evolution(...)

Run the time evolution for N_steps * dt.

TimeDependentTwoSiteTDVP.sweep([optimize])

One 'sweep' of a sweeper algorithm.

TimeDependentTwoSiteTDVP.switch_engine(...)

Initialize algorithm from another algorithm instance of a different class.

TimeDependentTwoSiteTDVP.update_env(...)

Do nothing; super().update_env() is called explicitly in update_local().

TimeDependentTwoSiteTDVP.update_local(theta, ...)

Perform algorithm-specific local update.

Class Attributes and Properties

TimeDependentTwoSiteTDVP.engine_params

TimeDependentTwoSiteTDVP.n_optimize

The number of sites to be optimized at once.

TimeDependentTwoSiteTDVP.time_dependent_H

whether the algorithm supports time-dependent H

TimeDependentTwoSiteTDVP.verbose

class tenpy.algorithms.tdvp.TimeDependentTwoSiteTDVP(psi, model, options, **kwargs)[source]

Bases: TimeDependentHAlgorithm, TwoSiteTDVPEngine

Variant of TwoSiteTDVPEngine that can handle time-dependent Hamiltonians.

See details in TimeDependentHAlgorithm as well.

reinit_model()[source]

Re-initialize a new model at current evolved_time.

Skips re-initialization if the model.options['time'] is the same as evolved_time. The model should read out the option 'time' and initialize the corresponding H(t).

EffectiveH[source]

alias of TwoSiteH

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

evolve(N_steps, dt)[source]

Evolve by N_steps * dt.

Parameters

N_steps (int) – The number of steps to evolve.

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 optimiztion 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 save psi, model and options 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 simulation from where we are now. It might contain an explicit copy of psi.

Return type

dict

get_sweep_schedule()[source]

Slightly different sweep schedule than DMRG

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. If None, 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 initializaiton.

  • 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 the MPSEnvironment; the psi 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 that resume_data['orthogonal_to'] takes precedence over the argument.

Options

Deprecated since version 0.6.0: Options LP, LP_age, RP and RP_age are now collected in a dictionary init_env_data with different keys init_LP, init_RP, age_LP, age_RP

Deprecated since version 0.8.0: Instead of passing the init_env_data as a option, it should be passed as dict entry of resume_data.

option Sweep.init_env_data: dict

Dictionary as returned by self.env.get_initialization_data() from get_initialization_data(). Deprecated, use the resume_data function/class argument instead.

option Sweep.orthogonal_to: list of MPS

Deprecated in favor of the orthogonal_to function argument (forwarded from the class argument) with the same effect.

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.

make_eff_H()[source]

Create new instance of self.EffectiveH at self.i0 and set it to self.eff_H.

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 the TwoSiteH and hence has n_optimize=2, while the SingleSiteDMRGEngine has n_optimize=1.

post_update_local(err, **update_data)[source]

Algorithm-specific actions to be taken after local update.

An example would be to collect statistics.

prepare_evolve(dt)[source]

Do nothing.

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

Array

reset_stats(resume_data=None)[source]

Reset the statistics. Useful if you want to start a new Sweep run.

This method is expected to be overwritten by subclass, and should then define self.update_stats and self.sweep_stats dicts consistent with the statistics generated by the algorithm particular to that subclass.

Parameters

resume_data (dict) – Given when resuming a simulation, as returned by get_resume_data(). Here, we read out the sweeps.

Options

option Sweep.chi_list: None | dict(int -> int)

By default (None) this feature is disabled. A dict allows to gradually increase the chi_max. An entry at_sweep: chi states that starting from sweep at_sweep, the value chi is to be used for trunc_params['chi_max']. For example chi_list={0: 50, 20: 100} uses chi_max=50 for the first 20 sweeps and chi_max=100 afterwards.

resume_run()[source]

Resume a run that was interrupted.

In case we saved an intermediate result at a checkpoint, this function allows to resume the run() 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 call run().

run()[source]

Perform a (real-)time evolution of psi by N_steps * dt.

You probably want to call this in a loop along with measurements. The recommended way to do this is via the RealTimeEvolution.

run_evolution(N_steps, dt)[source]

Run the time evolution for N_steps * dt.

Updates the model after each time step dt to account for changing H(t). For parameters see TimeEvolutionAlgorithm.

sweep(optimize=True)[source]

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

optimize (bool, optional) – Whether we actually optimize to find the ground state of the effective Hamiltonian. (If False, just update the environments).

Returns

max_trunc_err – Maximal truncation error introduced.

Return type

float

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 transfered. Another, but not too different algorithm subclass-class; e.g. you can switch from the TwoSiteDMRGEngine to the OneSiteDMRGEngine.

  • 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().

time_dependent_H = True

whether the algorithm supports time-dependent H

update_env(**update_data)[source]

Do nothing; super().update_env() is called explicitly in update_local().

update_local(theta, **kwargs)[source]

Perform algorithm-specific local update.

For two-site algorithms with n_optimize = 2, this always optimizes the sites i0 and i0 + 1. For single-site algorithms, the effective H only acts on site i0, but afterwards it also updates the bond to the right if move_right is True, or the bond to the left if move_right is False. Since the svd for truncation gives tensors to be multiplied into the tensors on both sides of the bond, tensors of two sites are updated even for single-site algorithms: when right-moving, site i0 + 1 is also updated; site i0 - 1 when left-moving.

Parameters

theta (Array) – Local single- or two-site wave function, as returned by prepare_update_local().

Returns

update_data – Data to be processed by update_env() and post_update_local(), e.g. containing the truncation error as err. If combine is set, it should also contain the U and VH from the SVD.

Return type

dict