SLOSBackend
The SLOSBackend
(for Strong Linear Optical Simulation) is a strong simulation backend that is able to compute efficiently
the entire output distribution by representing in memory a calculation path in which photons are added one by one,
with the best time complexity among strong simulation backends: \(\mathrm{O}(nC_n^{n+m-1})\). It is introduced in
Heurtel et al. [36].
The major downside of this backend is the memory intensive consumption, with the same complexity of \(\mathrm{O}(nC_n^{n+m-1})\). This backend is able to use masks to reduce the computation space, making it cheaper in memory and faster.
As such, this backend is well suited with a relatively small number of photons and modes (\(n, m < 20\)) when it is necessary to compute everything (or at least everything that befalls into a mask). If only a few output states are needed, other backends like NaiveBackend are more suited.
This backend is available in Processor by using the name "SLOS"
.
>>> import perceval as pcvl
>>> c = pcvl.Circuit(4) // pcvl.BS() // (2, pcvl.BS())
>>> backend = pcvl.SLOSBackend()
>>> backend.set_circuit(c)
>>> backend.set_input_state(pcvl.BasicState([1, 0, 1, 0]))
>>> print(backend.prob_distribution())
{
|1,0,1,0>: 0.2500000000000001
|1,0,0,1>: 0.2500000000000001
|0,1,1,0>: 0.2500000000000001
|0,1,0,1>: 0.2500000000000001
}
- class perceval.backends._slos.SLOSBackend(mask=None, use_symbolic=False)
- all_prob(input_state=None)
SLOS specific signature, to enhance optimization in some computations
- 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_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
- 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.