NaiveApproxBackend
Like the NaiveBackend, the NaiveApproxBackend
is a strong simulation backend that can compute a
single output probability amplitude at a time, by computing the permanent of a \(n \times n\) matrix.
However, instead of computing exactly the permanent, the NaiveApproxBackend
uses the Gurvit estimate algorithm
to approximate it with a 99% confidence interval (see Gurvits and Samorodnitsky [33]).
It is very efficient to compute very precise output states, but not to compute the whole distribution, and it can be used with more modes and photons than the NaiveBackend at the cost of losing precision on the result.
Thus, using it is not recommended in Simulator (except when using probability()
) or Processor,
but it is well suited for applications where only a few output probabilities matter with many photons.
If the whole or most of the computational space is needed, other backends like SLOSBackend are more suited.
This backend is available in Processor by using the name "NaiveApprox"
.
Unlike other backends, this backend needs a number of iterations to use to estimate the permanent. Also, in addition to the generic backend methods, this backend offers means to get a 99% confidence interval on the probability or a 99% sure error bound on the amplitude.
>>> import perceval as pcvl
>>> circuit_size = 60
>>> n_photons = 30
>>> backend = pcvl.NaiveApproxBackend(100_000_000) # Number of iterations; higher values reduce the error bound
>>> backend.set_circuit(pcvl.Unitary(pcvl.Matrix.random_unitary(circuit_size)))
>>> input_state = pcvl.BasicState([1]*n_photons + [0]*(circuit_size-n_photons))
>>> backend.set_input_state(input_state)
>>> interval = backend.probability_confidence_interval(BasicState([1]*n_photons + [0]*(circuit_size-n_photons)))
>>> print(f"Probability in {interval}")
Probability in [6.051670221391749e-20, 1.5297683283662674e-19]
- class perceval.backends._naive_approx.NaiveApproxBackend(gurvits_iterations=10000)
Naive algorithm with Gurvits computations of permanents
- Parameters:
gurvits_iterations (
int
) – Number of iterations to use for Gurvits estimation algorithm (default 10 000)
- all_prob(input_state=None)
Computes the list of probabilities of all states (respecting the mask if any was set). The order of the states can be retrieved using allstate_iterator()
- Return type:
list
[float
]
- clear_mask()
Removes any pre-existing mask
- evolve()
Evolves the input BasicState into a StateVector.
- Return type:
- property name: str
Returns the back-end name as a string
- prob_amplitude(output_state)
Computes the probability amplitude for a given output state. The input state and the circuit must already be set
- Return type:
complex
- prob_amplitude_with_error(output_state)
Computes the estimation of the probability amplitude along with an estimation of the 99% sure error bound.
- Return type:
tuple
[complex
,float
]
- prob_distribution()
Computes the probability distribution of all states (respecting the mask if any was set) under the form of a BSDistribution.
- Return type:
- probability(output_state)
Computes the probability for a given output state. The input state and the circuit must already be set
- Return type:
float
- probability_confidence_interval(output_state)
Computes the 99% confidence interval for the true value of the probability.
- Return type:
list
[float
]
- set_circuit(circuit)
Sets the circuit to simulate. This circuit must not contain polarized components (use PolarizationSimulator instead, if required).
- set_input_state(input_state)
Sets an input state for the simulation. This state has to be a Fock state without annotations.
- set_mask(masks, n=None)
Sets new masks, replacing the former ones if they exist. Clear possible cached data that depend on the mask. Masks are useful to limit strong simulation to only a part of the Fock space, ultimately saving memory and computation time.
- Parameters:
masks (str | list[str]) – Can be a mask or a list of masks. Each mask is expressed as a string where each character is a condition on one mode. Digits are fixing the number of photons whereas spaces or “*” are accepting any number of detections. e.g. using “****00” as a mask limits the simulation to output states ending in two empty modes.
n – The number of photons to instantiate the mask with. This corresponds to the total number of photons in your non-separated state.