netket.experimental.driver.TDVPSchmitt#
- class netket.experimental.driver.TDVPSchmitt[source]#
Bases:
TDVPBaseDriver
Variational time evolution based on the time-dependent variational principle which, when used with Monte Carlo sampling via
netket.vqs.MCState
, is the time-dependent VMC (t-VMC) method.This driver, which only works with standard MCState variational states, uses the regularization procedure described in M. Schmitt’s PRL 125 .
With the force vector
\[F_k=\langle \mathcal O_{\theta_k}^* E_{loc}^{\theta}\rangle_c\]and the quantum Fisher matrix
\[S_{k,k'} = \langle \mathcal O_{\theta_k} (\mathcal O_{\theta_{k'}})^*\rangle_c\]and for real parameters \(\theta\in\mathbb R\), the TDVP equation reads
\[q\big[S_{k,k'}\big]\theta_{k'} = -q\big[xF_k\big]\]Here, either \(q=\text{Re}\) or \(q=\text{Im}\) and \(x=1\) for ground state search or \(x=i\) (the imaginary unit) for real time dynamics.
For ground state search a regularization controlled by a parameter \(\rho\) can be included by increasing the diagonal entries and solving
\[q\big[(1+\rho\delta_{k,k'})S_{k,k'}\big]\theta_{k'} = -q\big[F_k\big]\]The TDVP class solves the TDVP equation by computing a pseudo-inverse of \(S\) via eigendecomposition yielding
\[S = V\Sigma V^\dagger\]with a diagonal matrix \(\Sigma_{kk}=\sigma_k\) Assuming that \(\sigma_1\) is the smallest eigenvalue, the pseudo-inverse is constructed from the regularized inverted eigenvalues
\[\tilde\sigma_k^{-1}=\frac{1}{\Big(1+\big(\frac{\epsilon_{SVD}}{\sigma_j/\sigma_1}\big)^6\Big)\Big(1+\big(\frac{\epsilon_{SNR}}{\text{SNR}(\rho_k)}\big)^6\Big)}\]with \(\text{SNR}(\rho_k)\) the signal-to-noise ratio of \(\rho_k=V_{k,k'}^{\dagger}F_{k'}\) (see [arXiv:1912.08828] for details).
- Inheritance
- __init__(operator, variational_state, ode_solver=None, *, t0=0.0, propagation_type='real', integrator=None, holomorphic=None, diag_shift=0.0, diag_scale=None, error_norm='qgt', rcond=1e-14, rcond_smooth=1e-08, snr_atol=1)[source]#
Initializes the time evolution driver.
- Parameters:
operator (
AbstractOperator
) – The generator of the dynamics (Hamiltonian for pure states, Lindbladian for density operators).variational_state (
VariationalState
) – The variational state.ode_solver (
AbstractSolver
) – Solving algorithm used the ODE.t0 (
float
) – Initial time at the start of the time evolution.propagation_type (
str
) – Determines the equation of motion: “real” for the real-time Schrödinger equation (SE), “imag” for the imaginary-time SE.error_norm (
str
|Callable
) – Norm function used to calculate the error with adaptive integrators. Can be either “euclidean” for the standard L2 vector norm \(w^\dagger w\), “maximum” for the maximum norm \(\max_i |w_i|\) or “qgt”, in which case the scalar product induced by the QGT \(S\) is used to compute the norm \(\Vert w \Vert^2_S = w^\dagger S w\) as suggested in PRL 125, 100503 (2020). Additionally, it possible to pass a custom function with signaturenorm(x: PyTree) -> float
which maps a PyTree of parametersx
to the corresponding norm. Note that norm is used in jax.jit-compiled code.holomorphic (
bool
|None
) – a flag to indicate that the wavefunction is holomorphic.diag_shift (
float
) – diagonal shift of the quantum geometric tensor (QGT)diag_scale (
float
|None
) – If not None rescales the diagonal shift of the QGTrcond (
float
) – Cut-off ratio for small singular \(\sigma_k\) values of the Quantum Geometric Tensor. For the purposes of rank determination, singular values are treated as zero if they are smaller than rcond times the largest singular valuesigma_{max}
.rcond_smooth (
float
) – Smooth cut-off ratio for singular values of the Quantum Geometric Tensor. This regularization parameter used with a similar effect to rcond but with a softer curve. See \(\epsilon_{SVD}\) in the formula above.snr_atol (
float
) – Noise regularisation absolute tolerance, meaning that eigenvalues of the S matrix that have a signal to noise ratio above this quantity will be (soft) truncated. This is \(\epsilon_{SNR}\) in the formulas above.integrator (AbstractSolver)
- Attributes
- dt#
Current time step.
- error_norm#
Returns the Callable function computing the error of the norm used for adaptive timestepping by the integrator.
Can be set to a Callable accepting a pytree and returning a real scalar, or a string between ‘euclidean’, ‘maximum’ or ‘qgt’.
- generator#
The generator of the dynamics as a function with signature generator(t: float) -> AbstractOperator
- integrator#
The underlying integrator which computes the time steps.
- ode_solver#
The underlying solver which solves the ODE at each time step.
- optimizer#
The optimizer used to update the parameters at every iteration.
- state#
Returns the machine that is optimized by this driver.
- step_count#
Returns a monotonic integer labelling all the steps performed by this driver. This can be used, for example, to identify the line in a log file.
- step_value#
- t#
Current time.
- t0#
The initial time set when the driver was created.
- Methods
- advance(T)[source]#
Advance the time propagation by
T
toself.t + T
.- Parameters:
T (
float
) – Length of the integration interval.
- estimate(observables)[source]#
Return MCMC statistics for the expectation value of observables in the current state of the driver.
- Parameters:
observables – A pytree of operators for which statistics should be computed.
- Returns:
A pytree of the same structure as the input, containing MCMC statistics for the corresponding operators as leaves.
- iter(T, *, tstops=None)[source]#
Returns a generator which advances the time evolution for an interval of length
T
, stopping attstops
.- Parameters:
T (
float
) – Length of the integration interval.tstops (
Sequence
[float
] |None
) – A sequence of stopping times, each within the interval[self.t0, self.t0 + T]
, at which this method will stop and yield. By default, a stop is performed after each time step (at potentially varying step size if an adaptive integrator is used).
- Yields:
The current step count.
- ode(t=None, w=None)[source]#
Evaluates the TDVP equation of motion
\[G(w) \dot w = \gamma F(w, t)\]where \(G(w)\) is the QGT, \(F(w, t)\) the gradient of
self.generator
and \(\gamma\) one of \(\gamma = -1\) (imaginary-time dynamics forMCState
), \(\gamma = -i\) (real-time dynamics forMCState
), or \(\gamma = 1\) (real-time dynamics forMCMixedState
).- Parameters:
t – Time (defaults to
self.t
).w – Variational parameters (defaults to
self.state.parameters
).
- Returns:
The time-derivative \(\dot w\).
- reset()[source]#
Resets the driver.
Subclasses should make sure to call
super().reset()
to ensure that the step count is set to 0.
- run(T, out=(), obs=None, *, tstops=None, show_progress=True, callback=None, timeit=False)[source]#
Runs the time evolution.
By default uses
netket.logging.JsonLog
. To know about the output format check it’s documentation. The logger object is also returned at the end of this function so that you can inspect the results without reading the json output.- Parameters:
T – The integration time period.
out – A logger object, or an iterable of loggers, to be used to store simulation log and data. If this argument is a string, it will be used as output prefix for the standard JSON logger.
obs – An iterable containing the observables that should be computed.
tstops – A sequence of stopping times, each within the interval
[self.t0, self.t0 + T]
, at which the driver will stop and perform estimation of observables, logging, and execute the callback function. By default, a stop is performed after each time step (at potentially varying step size if an adaptive integrator is used).show_progress – If true displays a progress bar (default=True)
callback – Callable or list of callable callback functions to be executed at each stopping time.
timeit (
bool
) – If True, provide timing information.