quimb.linalg.approx_spectral¶
Use stochastic Lanczos quadrature to approximate spectral function sums of any operator which has an efficient representation of action on a vector.
Attributes¶
Functions¶
|
A linear operator representing action of partially tracing a bipartite |
|
A linear operator representing action of partially tracing a tripartite |
|
Inner product between two vectors |
|
'Frobenius' norm of a vector. |
|
Calculate the approximate frobenius norm of any hermitian linear |
|
Generate a random array optionally orthogonal. |
|
Construct the tridiagonal lanczos matrix using only matvec operators. |
|
Find the eigen-values and -vectors of the Lanczos triadiagonal matrix. |
|
Spectral ritz function sum, weighted by ritz vectors. |
|
Extended percentile trimmed-mean. Makes the mean robust to asymmetric |
|
|
|
Simple standard deviation - don't invoke numpy for small lists. |
|
Make estimate by fitting exponential convergence to estimates. |
|
Make estimate from mean of last |
|
|
|
Get an estimate from samples. |
|
|
|
|
|
Approximate a spectral function, that is, the quantity |
|
|
|
|
|
|
|
|
|
|
|
Approximate the (Von Neumann) entropy of a pure state's subsystem. |
|
Approximate the trace sqrt of a pure state's subsystem. |
|
Estimate the norm of the partial transpose of a pure state's subsystem. |
|
Estimate the logarithmic negativity of a pure state's subsystem. |
|
Estimate the negativity of a pure state's subsystem. |
|
Generate a function that computes a spectral quantity of the subsystem |
Module Contents¶
- quimb.linalg.approx_spectral.reqs = '[cotengra,autoray]'¶
- quimb.linalg.approx_spectral.lazy_ptr_linop(psi_ab, dims, sysa, **linop_opts)[source]¶
A linear operator representing action of partially tracing a bipartite state, then multiplying another ‘unipartite’ state:
( | ) +-------+ | psi_a | ______ +_______+ / \ a| |b | +-------------+ | | psi_ab.H | | +_____________+ | | +-------------+ | | psi_ab | | +_____________+ | a| |b | | \______/
- quimb.linalg.approx_spectral.lazy_ptr_ppt_linop(psi_abc, dims, sysa, sysb, **linop_opts)[source]¶
A linear operator representing action of partially tracing a tripartite state, partially transposing the remaining bipartite state, then multiplying another bipartite state:
( | ) +--------------+ | psi_ab | +______________+ _____ a| ____ b| / \ | / a\ | |c | | | +-------------+ | | | | psi_abc.H | | \ / +-------------+ | X | / \ +-------------+ | | | | psi_abc | | | | +-------------+ | | \____/a |b |c | a| | \_____/
- Parameters:
psi_abc (ket) – State to partially trace, partially transpose, then dot with another ket, with size
prod(dims)
.prod(dims[sysa] + dims[sysb])
.dims (sequence of int) – The sub dimensions of
psi_abc
.sysa (int or sequence of int, optional) – Index(es) of the ‘a’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).sysa – Index(es) of the ‘b’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).
- quimb.linalg.approx_spectral.norm_fro_approx(A, **kwargs)[source]¶
Calculate the approximate frobenius norm of any hermitian linear operator:
\[\mathrm{Tr} \left[ A^{\dagger} A \right]\]- Parameters:
A (linear operator like) – Operator with a dot method, assumed to be hermitian, to estimate the frobenius norm of.
kwargs – Supplied to
approx_spectral_function()
.
- Return type:
- quimb.linalg.approx_spectral.random_rect(shape, dist='rademacher', orthog=False, norm=True, seed=False, dtype=complex)[source]¶
Generate a random array optionally orthogonal.
- quimb.linalg.approx_spectral.construct_lanczos_tridiag(A, K, v0=None, bsz=1, k_min=10, orthog=False, beta_tol=1e-06, seed=False, v0_opts=None)[source]¶
Construct the tridiagonal lanczos matrix using only matvec operators. This is a generator that iteratively yields the alpha and beta digaonals at each step.
- Parameters:
A (dense array, sparse matrix or linear operator) – The operator to approximate, must implement
.dot
method to compute its action on a vector.K (int, optional) – The maximum number of iterations and thus rank of the matrix to find.
v0 (vector, optional) – The starting vector to iterate with, default to random.
bsz (int, optional) – The block size (number of columns) of random vectors to iterate with.
k_min (int, optional) – The minimum size of the krylov subspace for form.
orthog (bool, optional) – If True, perform full re-orthogonalization for each new vector.
beta_tol (float, optional) – The ‘breakdown’ tolerance. If the next beta ceofficient in the lanczos matrix is less that this, implying that the full non-null space has been found, terminate early.
seed (bool, optional) – If True, seed the numpy random generator with a system random int.
- Yields:
alpha (sequence of float of length k) – The diagonal entries of the lanczos matrix.
beta (sequence of float of length k) – The off-diagonal entries of the lanczos matrix, with the last entry the ‘look’ forward value.
scaling (float) – How to scale the overall weights.
- quimb.linalg.approx_spectral.lanczos_tridiag_eig(alpha, beta, check_finite=True)[source]¶
Find the eigen-values and -vectors of the Lanczos triadiagonal matrix.
- quimb.linalg.approx_spectral.calc_trace_fn_tridiag(tl, tv, f, pos=True)[source]¶
Spectral ritz function sum, weighted by ritz vectors.
- quimb.linalg.approx_spectral.ext_per_trim(x, p=0.6, s=1.0)[source]¶
Extended percentile trimmed-mean. Makes the mean robust to asymmetric outliers, while using all data when it is nicely clustered. This can be visualized roughly as:
|--------|=========|--------| x x xx xx xxxxx xxx xx x x x
Where the inner range contains the central
p
proportion of the data, and the outer ranges entends this by a factor ofs
either side.- Parameters:
x (array) – Data to trim.
p (Proportion of data used to define the 'central' percentile.) – For example, p=0.5 gives the inter-quartile range.
s (Include data up to this factor times the central 'percentile' range) – away from the central percentile itself.
Returns
xt (array) – Trimmed data.
- quimb.linalg.approx_spectral.std(xs)[source]¶
Simple standard deviation - don’t invoke numpy for small lists.
- quimb.linalg.approx_spectral.calc_est_fit(estimates, conv_n, tau)[source]¶
Make estimate by fitting exponential convergence to estimates.
- quimb.linalg.approx_spectral.calc_est_window(estimates, conv_n)[source]¶
Make estimate from mean of last
m
samples, following:Take between
conv_n
and 12 estimates.Pair the estimates as they are alternate upper/lower bounds
Compute the standard error on the paired estimates.
- quimb.linalg.approx_spectral.single_random_estimate(A, K, bsz, beta_tol, v0, f, pos, tau, tol_scale, k_min=10, verbosity=0, *, seed=None, v0_opts=None, info=None, **lanczos_opts)[source]¶
- quimb.linalg.approx_spectral.calc_stats(samples, mean_p, mean_s, tol, tol_scale)[source]¶
Get an estimate from samples.
- quimb.linalg.approx_spectral.approx_spectral_function(A, f, tol=0.01, *, bsz=1, R=1024, R_min=3, tol_scale=1, tau=0.0001, k_min=10, k_max=512, beta_tol=1e-06, mpi=False, mean_p=0.7, mean_s=1.0, pos=False, v0=None, verbosity=0, single_precision='AUTO', info=None, progbar=False, plot=False, **lanczos_opts)[source]¶
Approximate a spectral function, that is, the quantity
Tr(f(A))
.- Parameters:
A (dense array, sparse matrix or LinearOperator) – Operator to approximate spectral function for. Should implement
A.dot(vec)
.f (callable) – Scalar function with which to act on approximate eigenvalues.
tol (float, optional) – Relative convergence tolerance threshold for error on mean of repeats. This can pretty much be relied on as the overall accuracy. See also
tol_scale
andtau
. Default: 1%.bsz (int, optional) – Number of simultenous vector columns to use at once, 1 equating to the standard lanczos method. If
bsz > 1
thenA
must implement matrix-matrix multiplication. This is a more performant way of essentially increasingR
, at the cost of more memory. Default: 1.R (int, optional) – The number of repeats with different initial random vectors to perform. Increasing this should increase accuracy as
sqrt(R)
. Cost of algorithm thus scales linearly withR
. Iftol
is non-zero, this is the maximum number of repeats.R_min (int, optional) – The minimum number of repeats to perform. Default: 3.
tau (float, optional) – The relative tolerance required for a single lanczos run to converge. This needs to be small enough that each estimate with a single random vector produces an unbiased sample of the operators spectrum..
k_min (int, optional) – The minimum size of the krylov subspace to form for each sample.
k_max (int, optional) – The maximum size of the kyrlov space to form. Cost of algorithm scales linearly with
K
. Iftau
is non-zero, this is the maximum size matrix to form.tol_scale (float, optional) – This sets the overall expected scale of each estimate, so that an absolute tolerance can be used for values near zero. Default: 1.
beta_tol (float, optional) – The ‘breakdown’ tolerance. If the next beta ceofficient in the lanczos matrix is less that this, implying that the full non-null space has been found, terminate early. Default: 1e-6.
mpi (bool, optional) – Whether to parallelize repeat runs over MPI processes.
mean_p (float, optional) – Factor for robustly finding mean and err of repeat estimates, see
ext_per_trim()
.mean_s (float, optional) – Factor for robustly finding mean and err of repeat estimates, see
ext_per_trim()
.v0 (vector, or callable) – Initial vector to iterate with, sets
R=1
if given. If callable, the function to produce a random intial vector (sequence).pos (bool, optional) – If True, make sure any approximate eigenvalues are positive by clipping below 0.
verbosity ({0, 1, 2}, optional) – How much information to print while computing.
single_precision ({'AUTO', False, True}, optional) – Try and convert the operator to single precision. This can lead to much faster operation, especially if a GPU is available. Additionally, double precision is not really needed given the stochastic nature of the algorithm.
lanczos_opts – Supplied to
single_random_estimate()
orconstruct_lanczos_tridiag()
.
- Returns:
The approximate value
Tr(f(a))
.- Return type:
scalar
See also
- quimb.linalg.approx_spectral.entropy_subsys_approx(psi_ab, dims, sysa, backend=None, **kwargs)[source]¶
Approximate the (Von Neumann) entropy of a pure state’s subsystem.
- Parameters:
psi_ab (ket) – Bipartite state to partially trace and find entopy of.
dims (sequence of int, optional) – The sub dimensions of
psi_ab
.sysa (int or sequence of int, optional) – Index(es) of the ‘a’ subsystem(s) to keep.
kwargs – Supplied to
approx_spectral_function()
.
- quimb.linalg.approx_spectral.tr_sqrt_subsys_approx(psi_ab, dims, sysa, backend=None, **kwargs)[source]¶
Approximate the trace sqrt of a pure state’s subsystem.
- Parameters:
psi_ab (ket) – Bipartite state to partially trace and find trace sqrt of.
dims (sequence of int, optional) – The sub dimensions of
psi_ab
.sysa (int or sequence of int, optional) – Index(es) of the ‘a’ subsystem(s) to keep.
kwargs – Supplied to
approx_spectral_function()
.
- quimb.linalg.approx_spectral.norm_ppt_subsys_approx(psi_abc, dims, sysa, sysb, backend=None, **kwargs)[source]¶
Estimate the norm of the partial transpose of a pure state’s subsystem.
- quimb.linalg.approx_spectral.logneg_subsys_approx(psi_abc, dims, sysa, sysb, **kwargs)[source]¶
Estimate the logarithmic negativity of a pure state’s subsystem.
- Parameters:
psi_abc (ket) – Pure tripartite state, for which estimate the entanglement between ‘a’ and ‘b’.
dims (sequence of int) – The sub dimensions of
psi_abc
.sysa (int or sequence of int, optional) – Index(es) of the ‘a’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).sysa – Index(es) of the ‘b’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).kwargs – Supplied to
approx_spectral_function()
.
- quimb.linalg.approx_spectral.negativity_subsys_approx(psi_abc, dims, sysa, sysb, **kwargs)[source]¶
Estimate the negativity of a pure state’s subsystem.
- Parameters:
psi_abc (ket) – Pure tripartite state, for which estimate the entanglement between ‘a’ and ‘b’.
dims (sequence of int) – The sub dimensions of
psi_abc
.sysa (int or sequence of int, optional) – Index(es) of the ‘a’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).sysa – Index(es) of the ‘b’ subsystem(s) to keep, with respect to all the dimensions,
dims
, (i.e. pre-partial trace).kwargs – Supplied to
approx_spectral_function()
.
- quimb.linalg.approx_spectral.gen_bipartite_spectral_fn(exact_fn, approx_fn, pure_default)[source]¶
Generate a function that computes a spectral quantity of the subsystem of a pure state. Automatically computes for the smaller subsystem, or switches to the approximate method for large subsystems.
- Parameters:
exact_fn (callable) – The function that computes the quantity on a density matrix, with signature:
exact_fn(rho_a, rank=...)
.approx_fn (callable) – The function that approximately computes the quantity using a lazy representation of the whole system. With signature
approx_fn(psi_ab, dims, sysa, **approx_opts)
.pure_default (float) – The default value when the whole state is the subsystem.
- Returns:
bipartite_spectral_fn – The function, with signature:
(psi_ab, dims, sysa, approx_thresh=2**13, **approx_opts)
- Return type:
callable