API
SLH rules
QuantumInputOutput.SLH — Type
SLHSLH triple with scattering matrix S, Lindblad term L and Hamiltonian H. S and L can also be vectors of scattering matrices and Linblad terms.
QuantumInputOutput.SLHqo — Type
SLHqoSLH triple where L and H are QuantumOptics.jl operators or time-dependent functions returning such operators. This is useful when you build models directly in a numeric basis, without symbolic translation.
QuantumInputOutput.get_scattering — Function
get_scattering(G)Return the scattering matrix S of an SLH or SLHqo object.
QuantumInputOutput.get_lindblad — Function
get_lindblad(G)Return the Lindblad vector L of an SLH or SLHqo object.
QuantumInputOutput.get_hamiltonian — Function
get_hamiltonian(G)Return the Hamiltonian H of an SLH or SLHqo object.
QuantumInputOutput.:▷ — Function
▷(G::SLH...)Creates a new SLH triple by cascading the SLH triples from first to last according to the rule
$SLH_1 \triangleright SLH_2 = ( S_2 S_1, L_2 + S_2 L_1, H_1 + H_2 - \frac{i}{2} L_2^\dagger S_2 L_1 - L_1^\dagger S_2^\dagger L_2 )$
Unicode \triangleright<tab> alias of cascade
▷(G::SLHqo...)Cascades SLHqo triples from first to last, allowing time-dependent L or H. If any input L or H is time-dependent, the result uses time-dependent functions.
QuantumInputOutput.cascade — Function
cascade(G::SLH...)Creates a new SLH triple by cascading the SLH triples from first to last according to the rule
\[SLH_1 \triangleright SLH_2 = ( S_2 S_1, L_2 + S_2 L_1, H_1 + H_2 - \frac{i}{2} L_2^\dagger S_2 L_1 - L_1^\dagger S_2^\dagger L_2 ) \]
See also ▷.
QuantumInputOutput.:⊞ — Function
⊞(G::SLH...)Creates a new SLH triple by concatenating the SLH triples according to the rule
\[SLH_1 \boxplus SLH_2 = \left( \begin{pmatrix} S_1 & 0 \\ 0 & S_2 \end{pmatrix}, \begin{pmatrix} L_1 \\ L_2 \end{pmatrix}, H_1 + H_2 \right)\]
Unicode \boxplus<tab> alias of concatenate
⊞(G::SLHqo...)Concatenates SLHqo triples, allowing time-dependent L or H. If any input L or H is time-dependent, the result uses time-dependent functions.
QuantumInputOutput.concatenate — Function
concatenate(G::SLH...)Creates a new SLH triple by concatenating the SLH triples according to the rule
\[SLH_1 \boxplus SLH_2 = \left( \begin{pmatrix} S_1 & 0 \\ 0 & S_2 \end{pmatrix}, \begin{pmatrix} L_1 \\ L_2 \end{pmatrix}, H_1 + H_2 \right)\]
See also ⊞.
QuantumInputOutput.feedback — Function
feedback(G::SLH, x::Int, y::Int)
feedback(G::SLH, connection::Pair{Int,Int})
feedback(G::SLH, connections::Pair{Int,Int}...)Apply the SLH feedback reduction rule to connect output port x to input port y. Multiple connections can be passed as x => y pairs; in that form the pairs are interpreted using the original port labels and reduced in sequence with the port bookkeeping handled automatically.
feedback(G::SLHqo, x::Int, y::Int)
feedback(G::SLHqo, connection::Pair{Int,Int})
feedback(G::SLHqo, connections::Pair{Int,Int}...)Apply the feedback reduction rule to an SLHqo triple. Time-dependent L or H are preserved in the reduced model.
Translation
QuantumInputOutput.translate_qo — Function
translate_qo(op, b::QuantumOpticsBase.Basis; parameter=Dict(), time_parameter=Dict(),
level_map=nothing, operators=Dict(), op_type=sparse)Translate a symbolic operator op into a numeric QuantumOptics.jl operator with the corresponding basis b. The dictionary parameter substitutes symbolic parameters with numbers. Time-dependent functions can be provide with the dictionary time_parameter. If time_parameter is non-empty, the result is a time-dependent function t -> op(t). The kwarg level_map=nothing is used to provide the names of levels for transition operators. The operator type which should be returned can be set with the kwarg op_type=sparse and a list of user-defined operators (e.g. on a different basis than b) can be provided with the dictionary operators=Dict(). These operators will then be used to replace the symbolic expressions.
Pulses
QuantumInputOutput.u_to_gu — Function
u_to_gu(u, T)Compute the virtual-cavity coupling $g_u(t)$ from an input mode u(t) sampled on T with a linear interpolation. Returns a callable t -> g_u(t).
QuantumInputOutput.v_to_gv — Function
v_to_gv(v, T)Compute the virtual-cavity coupling $g_v(t)$ from an output mode v(t) sampled on T with a linear interpolation. Returns a callable t -> g_v(t).
QuantumInputOutput.u_to_gu_Gauss — Function
u_to_gu_Gauss(τ, σ; δ=0)Compute $g_u(t)$ for a Gaussian input mode u(t) with delay τ, width σ, and detuning δ:
$u(t) = 1/\sqrt{ \sigma \sqrt{\pi} } e^{-(t - \tau)^2 / 2 \sigma^2 )} e^{i \delta t}$
Returns a callable t -> g_u(t).
QuantumInputOutput.v_to_gv_Gauss — Function
v_to_gv_Gauss(τ, σ; δ=0)Compute $g_v(t)$ for a Gaussian output mode v(t) with delay τ, width σ, and detuning δ:
$v(t) = 1/\sqrt{ \sigma \sqrt{\pi} } e^{-(t - \tau)^2 / 2 \sigma^2 )} e^{i \delta t}$
Returns a callable t -> g_v(t).
QuantumInputOutput.u_eff — Function
u_eff(u_fcts, gu_fcts, T, i)
u_eff(u_fcts, T, i)Compute the effective input mode $u_i^{\mathrm{eff}}(t)$ for a system with multiple input modes, due to the pulse distortion from the subsequent input cavities. The returned mode function is the relevant one for $g_{u_i}$. The input modes in the list u_fcts need to be sorted starting with the first input cavity before the system.
All kwargs are passed on to the ODE solver.
QuantumInputOutput.v_eff — Function
v_eff(v_fcts, gv_fcts, T, i)
v_eff(v_fcts, T, i)Compute the effective output mode $v_i^{\mathrm{eff}}(t)$ for a system with multiple output modes, due to the pulse distortion from the preceding output cavities. The returned mode function is the relevant one for $g_{v_i}$. The output modes in the list v_fcts need to be sorted starting with the first output cavity after the system.
All kwargs are passed on to the ODE solver.
QuantumInputOutput.uv_to_gin — Function
uv_to_gin(u, v, T)Compute the in-coupling strength $\tilde g_{\mathrm{in},v,u}(t)$ for a delay cavity that absorbs an incoming pulse v(t) while simultaneously emitting the desired pulse u(t). Returns callable t -> g_in(t) based on samples on T with a linear interpolation.
QuantumInputOutput.uv_to_gout — Function
uv_to_gout(u, v, T)Compute the out-coupling strength $\tilde g_{\mathrm{out},u,v}(t)$ for a delay cavity that absorbs an incoming pulse v(t) while simultaneously emitting the desired pulse u(t). Returns callable t -> g_out(t) based on samples on T with a linear interpolation.
Interaction Picture
QuantumInputOutput.interaction_picture_M — Function
interaction_picture_M(A, T; alg = Tsit5(), kwargs...)Solve the Heisenberg coefficient-matrix equation for an interaction picture, dM/dt = A(t) * M(t) with M(0) = I, and return the ODE solution.
The function A(t) must return an N × N complex matrix describing the linear mode coupling in the interaction picture. The returned solution supports interpolation M(t) for continuous t.
QuantumInputOutput.interaction_picture_M_2modes_equal — Function
interaction_picture_M_2modes_equal(u, T)Analytic interaction-picture coefficient matrix for two modes when u(t) = v(t). The matrix is
\[M(t) = \begin{bmatrix} \cos \theta(t) & -\sin \theta(t) \\ \sin \theta(t) & \cos \theta(t) \end{bmatrix},\]
where
\[\sin^2 \theta(t) = \int_0^t |u(t')|^2\,dt'.\]
QuantumInputOutput.interaction_picture_A_2modes — Function
interaction_picture_A_2modes(g1, g2)Coefficient matrix A(t) for two virtual modes (1, 2),
\[A(t) = \frac{1}{2} \begin{bmatrix} 0 & g_1(t) g_2^*(t) \\ -g_1^*(t) g_2(t) & 0 \end{bmatrix}\]
All couplings may be time-dependent or constant.
QuantumInputOutput.interaction_picture_A_3modes — Function
interaction_picture_A_3modes(g1, g2, g3)Coefficient matrix A(t) for three interacting modes ordered as (1, 2, 3)
\[A(t) = \frac{1}{2} \begin{bmatrix} 0 & g_1(t) g_2^*(t) & g_1(t) g_3^*(t) \\ -g_1^*(t) g_2(t) & 0 & g_2(t) g_3^*(t) \\ -g_1^*(t) g_3(t) & -g_2^*(t) g_3(t) & 0 \end{bmatrix}\]
All couplings may be time-dependent or constant.
QuantumInputOutput.interaction_picture_A_4modes — Function
interaction_picture_A_4modes(g1, g2, g3, g4)Coefficient matrix A(t) for four interacting modes ordered as (1, 2, 3, 4)
\[A(t) = \frac{1}{2} \begin{bmatrix} 0 & g_1(t) g_2^*(t) & g_1(t) g_3^*(t) & g_1(t) g_4^*(t) \\ -g_1^*(t) g_2(t) & 0 & g_2(t) g_3^*(t) & g_2(t) g_4^*(t) \\ -g_1^*(t) g_3(t) & -g_2^*(t) g_3(t) & 0 & g_3(t) g_4^*(t) \\ -g_1^*(t) g_4(t) & -g_2^*(t) g_4(t) & -g_3^*(t) g_4(t) & 0 \end{bmatrix}\]
All couplings may be time-dependent or constant.
Correlations
QuantumInputOutput.two_time_corr_matrix — Function
two_time_corr_matrix(T, ρt, f, Ls; kwargs...)
two_time_corr_matrix(T, ρt, H, J, Ls; kwargs...)Compute the two-time correlation matrix $g^{(1)}(t_1, t_2) = \langle L_s^\dagger(t_1) L_s(t_2) \rangle$ on the time grid T for the operator Ls. The first method supports time-dependent generators with the function f; the second is for time-independent H and J.
Utilities
QuantumInputOutput.substitute_operators — Function
substitute_operators(op, dict::Dict; replace_adjoint=true)Like substitute(op, dict::Dict) but with special handling for QMul and QAdd. This is needed if an operator is substitute by a QMul or QAdd, e.g. $a_1 \rightarrow g_2 a_2 + g_3 a_3$
If replace_adjoint=true, the dictionary is extended with adjoint substitutions for all key/value pairs, i.e. adjoint(key) => adjoint(value).