tensorcircuit.timeevol¶
Analog time evolution engines
- tensorcircuit.timeevol.chebyshev_evol(hamiltonian: Any, initial_state: Any, t: float, spectral_bounds: Tuple[float, float], k: int, M: int) Any[source]¶
Chebyshev evolution method by expanding the time evolution exponential operator in Chebyshev series. Note the state returned is not normalized. But the norm should be very close to 1 for sufficiently large k and M, which can serve as a accuracy check of the final result.
- Parameters:
hamiltonian (Any) – Hamiltonian matrix (sparse or dense)
initial_state (Tensor) – Initial state vector
time (float) – Time to evolve
spectral_bounds (Tuple[float, float]) – Spectral bounds for the Hamiltonian (Emax, Emin)
k (int) – Number of Chebyshev coefficients, a good estimate is k > t*(Emax-Emin)/2
M (int) – Number of iterations to estimate Bessel function, a good estimate is given by estimate_M helper method.
- Returns:
Evolved state
- Return type:
Tensor
- tensorcircuit.timeevol.ed_evol(h: Any, psi0: Any, tlist: Any, callback: Callable[[...], Any] | None = None) Any¶
Fast implementation of time independent Hamiltonian evolution using eigendecomposition. By default, performs imaginary time evolution.
- Parameters:
h (Tensor) – Time-independent Hamiltonian matrix
hamiltonian (Tensor) – alias for the argument
hpsi0 (Tensor) – Initial state vector
initial_state (Tensor) – alias for the argument
psi0tlist (Tensor) – Time points for evolution
times (Tensor) – alias for the argument
tlistcallback (Optional[Callable[..., Any]], optional) – Optional function to process state at each time point
- Returns:
Evolution results at each time point. If callback is None, returns state vectors; otherwise returns callback results
- Return type:
Tensor
- Example:
>>> import tensorcircuit as tc >>> import numpy as np >>> # Define a simple 2-qubit Hamiltonian >>> h = tc.array_to_tensor([ ... [1.0, 0.0, 0.0, 0.0], ... [0.0, -1.0, 2.0, 0.0], ... [0.0, 2.0, -1.0, 0.0], ... [0.0, 0.0, 0.0, 1.0] ... ]) >>> # Initial state |00> >>> psi0 = tc.array_to_tensor([1.0, 0.0, 0.0, 0.0]) >>> # Evolution times >>> times = tc.array_to_tensor([0.0, 0.5, 1.0]) >>> # Evolve and get states >>> states = tc.experimental.hamiltonian_evol(times, h, psi0) >>> print(states.shape) # (3, 4)
Note
The Hamiltonian must be time-independent
For time-dependent Hamiltonians, use
evol_localorevol_globalinsteadThe evolution is performed in imaginary time by default (factor -t in exponential)
The state is automatically normalized at each time point
- tensorcircuit.timeevol.estimate_M(t: float, spectral_bounds: Tuple[float, float], k: int) int[source]¶
estimate M for Bessel function iterations
- Parameters:
t (float) – time
spectral_bounds (Tuple[float, float]) – spectral bounds (Emax, Emin)
k (int) – k
- Returns:
M
- Return type:
int
- tensorcircuit.timeevol.estimate_k(t: float, spectral_bounds: Tuple[float, float]) int[source]¶
estimate k for chebyshev expansion
- Parameters:
t (float) – time
spectral_bounds (Tuple[float, float]) – spectral bounds (Emax, Emin)
- Returns:
k
- Return type:
int
- tensorcircuit.timeevol.estimate_spectral_bounds(h: Any, n_iter: int = 30, psi0: Any | None = None) Tuple[float, float][source]¶
Lanczos algorithm to estimate the spectral bounds of a Hamiltonian. Just for quick run before chebyshev_evol, non jit-able.
- Parameters:
h (Any) – Hamiltonian matrix.
n_iter (int) – iteration number.
psi0 (Optional[Any]) – Optional initial state.
- Returns:
(E_max, E_min)
- tensorcircuit.timeevol.evol_global(c: Any, h_fun: Callable[[...], Any], t: float, *args: Any, **solver_kws: Any) Any[source]¶
ode evolution of time dependent Hamiltonian on circuit of all qubits [only jax backend support for now]
- Parameters:
c (Circuit) – _description_
h_fun (Callable[..., Tensor]) – h_fun should return a SPARSE Hamiltonian matrix with input arguments
timeand*argshamiltonian (Callable[..., Tensor]) – alias for the argument
h_funt (float) – _description_
times (float) – alias for the argument
t
- Returns:
_description_
- Return type:
- tensorcircuit.timeevol.evol_local(c: Any, index: Sequence[int], h_fun: Callable[[...], Any], t: float, *args: Any, **solver_kws: Any) Any[source]¶
ode evolution of time dependent Hamiltonian on circuit of given indices [only jax backend support for now]
- Parameters:
c (Circuit) – _description_
index (Sequence[int]) – qubit sites to evolve
h_fun (Callable[..., Tensor]) – h_fun should return a dense Hamiltonian matrix with input arguments
timeand*argshamiltonian (Callable[..., Tensor]) – alias for the argument
h_funt (float) – evolution time
times (float) – alias for the argument
t
- Returns:
_description_
- Return type:
- tensorcircuit.timeevol.hamiltonian_evol(h: Any, psi0: Any, tlist: Any, callback: Callable[[...], Any] | None = None) Any[source]¶
Fast implementation of time independent Hamiltonian evolution using eigendecomposition. By default, performs imaginary time evolution.
- Parameters:
h (Tensor) – Time-independent Hamiltonian matrix
hamiltonian (Tensor) – alias for the argument
hpsi0 (Tensor) – Initial state vector
initial_state (Tensor) – alias for the argument
psi0tlist (Tensor) – Time points for evolution
times (Tensor) – alias for the argument
tlistcallback (Optional[Callable[..., Any]], optional) – Optional function to process state at each time point
- Returns:
Evolution results at each time point. If callback is None, returns state vectors; otherwise returns callback results
- Return type:
Tensor
- Example:
>>> import tensorcircuit as tc >>> import numpy as np >>> # Define a simple 2-qubit Hamiltonian >>> h = tc.array_to_tensor([ ... [1.0, 0.0, 0.0, 0.0], ... [0.0, -1.0, 2.0, 0.0], ... [0.0, 2.0, -1.0, 0.0], ... [0.0, 0.0, 0.0, 1.0] ... ]) >>> # Initial state |00> >>> psi0 = tc.array_to_tensor([1.0, 0.0, 0.0, 0.0]) >>> # Evolution times >>> times = tc.array_to_tensor([0.0, 0.5, 1.0]) >>> # Evolve and get states >>> states = tc.experimental.hamiltonian_evol(times, h, psi0) >>> print(states.shape) # (3, 4)
Note
The Hamiltonian must be time-independent
For time-dependent Hamiltonians, use
evol_localorevol_globalinsteadThe evolution is performed in imaginary time by default (factor -t in exponential)
The state is automatically normalized at each time point
- tensorcircuit.timeevol.krylov_evol(hamiltonian: Any, initial_state: Any, times: Any, subspace_dimension: int, callback: Callable[[Any], Any] | None = None, scan_impl: bool = False) Any[source]¶
Perform quantum state time evolution using Krylov subspace method.
- Parameters:
hamiltonian (Tensor) – Sparse or dense Hamiltonian matrix
initial_state (Tensor) – Initial quantum state
times (Tensor) – List of time points
subspace_dimension (int) – Krylov subspace dimension
callback (Optional[Callable[[Any], Any]], optional) – Optional callback function applied to quantum state at each evolution time point, return some observables
scan_impl (bool, optional) – whether use scan implementation, suitable for jit but may be slow on numpy defaults False, True not work for tensorflow backend + jit, due to stupid issue of tensorflow context separation and the notorious inaccesibletensor error
- Returns:
List of evolved quantum states, or list of callback function results (if callback provided)
- Return type:
Any
- tensorcircuit.timeevol.lanczos_iteration(hamiltonian: Any, initial_vector: Any, subspace_dimension: int) Tuple[Any, Any][source]¶
Use Lanczos algorithm to construct orthogonal basis and projected Hamiltonian of Krylov subspace.
- Parameters:
hamiltonian (Tensor) – Sparse or dense Hamiltonian matrix
initial_vector (Tensor) – Initial quantum state vector
subspace_dimension (int) – Dimension of Krylov subspace
- Returns:
Tuple containing (basis matrix, projected Hamiltonian)
- Return type:
Tuple[Tensor, Tensor]
- tensorcircuit.timeevol.lanczos_iteration_scan(hamiltonian: Any, initial_vector: Any, subspace_dimension: int) Tuple[Any, Any][source]¶
Use Lanczos algorithm to construct orthogonal basis and projected Hamiltonian of Krylov subspace, using tc.backend.scan for JIT compatibility.
- Parameters:
hamiltonian (Tensor) – Sparse or dense Hamiltonian matrix
initial_vector (Tensor) – Initial quantum state vector
subspace_dimension (int) – Dimension of Krylov subspace
- Returns:
Tuple containing (basis matrix, projected Hamiltonian)
- Return type:
Tuple[Tensor, Tensor]
- tensorcircuit.timeevol.ode_evol_global(hamiltonian: Callable[[...], Any], initial_state: Any, times: Any, callback: Callable[[...], Any] | None = None, *args: Any, **solver_kws: Any) Any[source]¶
ODE-based time evolution for a time-dependent Hamiltonian acting on the entire system. This function solves the time-dependent Schrodinger equation using numerical ODE integration. The Hamiltonian is applied to the full system and should be provided in sparse matrix format for efficiency. The ode_backend parameter defaults to ‘jaxode’ (which uses
jax.experimental.ode.odeintwith a default solver of ‘Dopri5’); if set to ‘diffrax’, it usesdiffrax.diffeqsolveinstead (with a default solver of ‘Tsit5’).Note: This function currently only supports the JAX backend.
- Parameters:
hamiltonian (Callable[..., Tensor]) – A function that returns a sparse Hamiltonian matrix for the full system. The function signature should be
hamiltonian(time, *args) -> Tensor.initial_state (Tensor) – The initial quantum state vector.
times (Tensor) – Time points for which to compute the evolution. Should be a 1D array of times.
callback (Optional[Callable[..., Tensor]]) – Optional function to apply to the state at each time step.
args (tuple | list) – Additional arguments to pass to the Hamiltonian function.
solver_kws (dict) –
Additional keyword arguments to pass to the ODE solver.
ode_backend='jaxode'(default) usesjax.experimental.ode.odeint;ode_backend='diffrax'usesdiffrax.diffeqsolve.rtol(default: 1e-8) andatol(default: 1e-8) are used to determine how accurately you would like the numerical approximation to your equation.The
solverparameter accepts one of {‘Tsit5’ (default), ‘Dopri5’, ‘Dopri8’, ‘Kvaerno5’} and only works whenode_backend='diffrax'.t0(default: 0.01) specifies the initial step size and only works whenode_backend='diffrax'.max_steps(default: 4096) The maximum number of steps to take before quitting the computation unconditionally and only works whenode_backend='diffrax'.
- Returns:
Evolved quantum states at the specified time points. If callback is provided, returns the callback results; otherwise returns the state vectors.
- Return type:
Tensor
- tensorcircuit.timeevol.ode_evol_local(hamiltonian: Callable[[...], Any], initial_state: Any, times: Any, index: Sequence[int], callback: Callable[[...], Any] | None = None, *args: Any, **solver_kws: Any) Any[source]¶
ODE-based time evolution for a time-dependent Hamiltonian acting on a subsystem of qubits. This function solves the time-dependent Schrodinger equation using numerical ODE integration. The Hamiltonian is applied only to a specific subset of qubits (indices) in the system. The ode_backend parameter defaults to ‘jaxode’ (which uses
jax.experimental.ode.odeintwith a default solver of ‘Dopri5’); if set to ‘diffrax’, it usesdiffrax.diffeqsolveinstead (with a default solver of ‘Tsit5’).Note: This function currently only supports the JAX backend.
- Parameters:
hamiltonian (Callable[..., Tensor]) – A function that returns a dense Hamiltonian matrix for the specified subsystem size. The function signature should be
hamiltonian(time, *args) -> Tensor.initial_state (Tensor) – The initial quantum state vector of the full system.
times (Tensor) – Time points for which to compute the evolution. Should be a 1D array of times.
index (Sequence[int]) – Indices of qubits where the Hamiltonian is applied.
callback (Optional[Callable[..., Tensor]]) – Optional function to apply to the state at each time step.
args – Additional arguments to pass to the Hamiltonian function.
solver_kws (dict) –
Additional keyword arguments to pass to the ODE solver.
ode_backend='jaxode'(default) usesjax.experimental.ode.odeint;ode_backend='diffrax'usesdiffrax.diffeqsolve.rtol(default: 1e-8) andatol(default: 1e-8) are used to determine how accurately you would like the numerical approximation to your equation.The
solverparameter accepts one of {‘Tsit5’ (default), ‘Dopri5’, ‘Dopri8’, ‘Kvaerno5’} and only works whenode_backend='diffrax'.t0(default: 0.01) specifies the initial step size and only works whenode_backend='diffrax'.max_steps(default: 4096) The maximum number of steps to take before quitting the computation unconditionally and only works whenode_backend='diffrax'.
- Returns:
Evolved quantum states at the specified time points. If callback is provided, returns the callback results; otherwise returns the state vectors.
- Return type:
Tensor